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.
Changeset 2715 for trunk/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modbcfunction.F – NEMO

Ignore:
Timestamp:
2011-03-30T17:58:35+02:00 (13 years ago)
Author:
rblod
Message:

First attempt to put dynamic allocation on the trunk

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modbcfunction.F

    r2528 r2715  
    3535      Use Agrif_Update 
    3636      Use Agrif_fluxmod 
     37      Use Agrif_Save 
    3738C              
    3839      IMPLICIT NONE 
     
    6566     &                     Agrif_Init_variable1d, 
    6667     &                     Agrif_Init_variable2d, 
    67      &                     Agrif_Init_variable3d 
     68     &                     Agrif_Init_variable3d, 
     69     &                     Agrif_Init_variable4d 
    6870      end interface        
    6971C 
     
    7678     &                     Agrif_update_var5d 
    7779      end interface        
     80       
     81      interface Agrif_Save_Forrestore 
     82         module procedure Agrif_Save_Forrestore0d,     
     83     &                    Agrif_Save_Forrestore2d,     
     84     &                    Agrif_Save_Forrestore3d,     
     85     &                    Agrif_Save_Forrestore4d 
     86      end interface 
    7887C 
    7988      Contains  
     
    255264      LOGICAL, OPTIONAL :: Interpolationshouldbemade 
    256265C 
    257       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     266      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    258267      TYPE(Agrif_PVariable),Pointer ::tabvars 
    259268       
     
    265274C      
    266275 
    267       if (tabvarsindic <=0) then 
    268       tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) 
    269       else 
    270       tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) 
     276      indic = tabvarsindic 
     277      if (tabvarsindic >=0) then 
     278        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     279          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     280        endif 
     281      endif 
     282       
     283      if (indic <=0) then 
     284      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     285      else 
     286      tabvars=>Agrif_Curgrid % tabvars(indic) 
    271287      endif   
    272288       
     
    307323      INTEGER, OPTIONAL      :: interp,interp1,interp2,interp3 
    308324C 
    309       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     325      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
     326      TYPE(Agrif_PVariable),Pointer ::tabvars 
     327       
     328     
     329C 
    310330C 
    311331C     Begin  
    312332C 
    313       Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp =  
     333C      
     334      indic = tabvarsindic 
     335      if (tabvarsindic >=0) then 
     336        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     337          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     338        endif 
     339      endif  
     340       
     341      if (indic <=0) then 
     342      tabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     343      else 
     344      tabvars=>Agrif_Mygrid % tabvars(indic) 
     345      endif       
     346C 
     347C     Begin  
     348C 
     349      tabvars % var % Typeinterp =  
    314350     &    Agrif_Constant 
    315351      IF (present(interp)) THEN 
    316       Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp =  
     352      tabvars % var % Typeinterp =  
    317353     &           interp 
    318354      ENDIF 
    319355      IF (present(interp1)) THEN 
    320       Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp(1) =  
     356      tabvars % var % Typeinterp(1) =  
    321357     &           interp1 
    322358      ENDIF 
    323359      IF (present(interp2)) THEN 
    324       Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp(2) =  
     360      tabvars % var % Typeinterp(2) =  
    325361     &           interp2 
    326362      ENDIF 
    327363      IF (present(interp3)) THEN 
    328       Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp(3) =  
     364      tabvars % var % Typeinterp(3) =  
    329365     &           interp3 
    330366      ENDIF 
     
    353389      INTEGER, OPTIONAL      :: interp11,interp12,interp21,interp22 
    354390C 
    355       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     391      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    356392      TYPE(Agrif_PVariable),Pointer ::tabvars 
    357393       
     
    363399C      
    364400 
    365       if (tabvarsindic <=0) then 
    366       tabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) 
    367       else 
    368       tabvars=>Agrif_Mygrid % tabvars(tabvarsindic) 
     401      indic = tabvarsindic 
     402      if (tabvarsindic >=0) then 
     403        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     404          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     405        endif 
     406      endif  
     407       
     408      if (indic <=0) then 
     409      tabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     410      else 
     411      tabvars=>Agrif_Mygrid % tabvars(indic) 
    369412      endif  
    370413C 
     
    460503     &       update2, update3,update4,update5 
    461504C 
    462       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     505      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
     506      TYPE(Agrif_PVariable),Pointer :: roottabvars       
    463507C 
    464508C 
    465509C     Begin  
    466 C 
    467       Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate =  
     510 
     511      indic = tabvarsindic 
     512 
     513      if (tabvarsindic >=0) then 
     514        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     515          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     516        endif 
     517      endif  
     518       
     519      if (indic <=0) then 
     520      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic)       
     521      else 
     522      roottabvars => Agrif_Mygrid % tabvars(indic) 
     523      endif 
     524       
     525C 
     526      roottabvars% var % typeupdate =  
    468527     &                   Agrif_Update_Copy 
    469528       
    470529      IF (present(update)) THEN 
    471         Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate =  
     530        roottabvars% var % typeupdate =  
    472531     &           update 
    473532      ENDIF 
    474533      IF (present(update1)) THEN 
    475         Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(1) =  
     534        roottabvars% var % typeupdate(1) =  
    476535     &           update1 
    477536      ENDIF   
    478537      IF (present(update2)) THEN 
    479         Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(2) =  
     538        roottabvars% var % typeupdate(2) =  
    480539     &           update2 
    481540      ENDIF   
    482541      IF (present(update3)) THEN 
    483         Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(3) =  
     542        roottabvars% var % typeupdate(3) =  
    484543     &           update3 
    485544      ENDIF  
    486545      IF (present(update4)) THEN 
    487         Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(4) =  
     546        roottabvars% var % typeupdate(4) =  
    488547     &           update4 
    489548      ENDIF        
    490549      IF (present(update5)) THEN 
    491         Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(5) =  
     550        roottabvars% var % typeupdate(5) =  
    492551     &           update5 
    493552      ENDIF                   
     
    513572C     Arguments       
    514573C 
    515       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     574      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
    516575C 
    517576C     Begin  
    518577C 
    519 C 
    520       Agrif_Mygrid%tabvars(tabvarsindic)%var % restaure = .TRUE. 
     578      indic = tabvarsindic 
     579      if (tabvarsindic >=0) then 
     580        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     581          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     582        endif 
     583      endif   
     584C 
     585      Agrif_Mygrid%tabvars(indic)%var % restaure = .TRUE. 
    521586C 
    522587      End Subroutine Agrif_Set_restore 
     
    530595 
    531596      INTEGER :: tabvarsindic0 ! indice of the variable in tabvars 
    532       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     597      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
    533598      External :: procname 
    534599      Optional ::  procname 
     
    536601      if (Agrif_Root()) Return 
    537602C       
     603      indic = tabvarsindic 
     604      if (tabvarsindic >=0) then 
     605        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     606          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     607        endif 
     608      endif  
     609       
    538610      if (present(procname)) then 
    539       CALL Agrif_Interp_variable(tabvarsindic0,tabvarsindic,procname) 
    540       CALL Agrif_Bc_variable(tabvarsindic0,tabvarsindic,1.,procname) 
    541       else 
    542       CALL Agrif_Interp_variable(tabvarsindic0,tabvarsindic) 
    543       CALL Agrif_Bc_variable(tabvarsindic0,tabvarsindic,1.) 
     611      CALL Agrif_Interp_variable(tabvarsindic0,indic,procname) 
     612      CALL Agrif_Bc_variable(tabvarsindic0,indic,1.,procname) 
     613      else 
     614      CALL Agrif_Interp_variable(tabvarsindic0,indic) 
     615      CALL Agrif_Bc_variable(tabvarsindic0,indic,1.) 
    544616      endif 
    545617 
     
    553625 
    554626      REAL, DIMENSION(:) :: q 
    555       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     627      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    556628      External :: procname 
    557629      Optional ::  procname 
     
    559631C 
    560632      if (Agrif_Root()) Return 
     633C       
     634      indic = tabvarsindic 
     635      if (tabvarsindic >=0) then 
     636        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     637          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     638        endif 
     639      endif       
    561640C 
    562641      if (present(procname)) then 
    563       CALL Agrif_Interp_variable(q,tabvarsindic,procname) 
    564       CALL Agrif_Bc_variable(q,tabvarsindic,1.,procname) 
    565       else 
    566       CALL Agrif_Interp_variable(q,tabvarsindic) 
    567       CALL Agrif_Bc_variable(q,tabvarsindic,1.) 
     642      CALL Agrif_Interp_variable(q,indic,procname) 
     643      CALL Agrif_Bc_variable(q,indic,1.,procname) 
     644      else 
     645      CALL Agrif_Interp_variable(q,indic) 
     646      CALL Agrif_Bc_variable(q,indic,1.) 
    568647      endif 
    569648 
     
    579658      External :: procname 
    580659      Optional ::  procname 
     660      integer :: indic 
    581661 
    582662C 
    583663      if (Agrif_Root()) Return 
    584664C 
     665      indic = tabvarsindic 
     666      if (tabvarsindic >=0) then 
     667        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     668          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     669        endif 
     670      endif  
     671       
    585672      if (present(procname)) then 
    586       CALL Agrif_Interp_variable(q,tabvarsindic,procname) 
    587       CALL Agrif_Bc_variable(q,tabvarsindic,1.,procname) 
    588       else 
    589       CALL Agrif_Interp_variable(q,tabvarsindic) 
    590       CALL Agrif_Bc_variable(q,tabvarsindic,1.) 
     673      CALL Agrif_Interp_variable(q,indic,procname) 
     674      CALL Agrif_Bc_variable(q,indic,1.,procname) 
     675      else 
     676      CALL Agrif_Interp_variable(q,indic) 
     677      CALL Agrif_Bc_variable(q,indic,1.) 
    591678      endif 
    592679 
     
    601688 
    602689      REAL,  DIMENSION(:,:,:) :: q 
    603       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     690      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
    604691      External :: procname 
    605692      Optional ::  procname 
    606693C 
    607694      if (Agrif_Root()) Return 
     695C       
     696      indic = tabvarsindic 
     697      if (tabvarsindic >=0) then 
     698        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     699          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     700        endif 
     701      endif       
    608702C 
    609703      if (present(procname)) then 
    610       CALL Agrif_Interp_variable(q,tabvarsindic,procname) 
    611       CALL Agrif_Bc_variable(q,tabvarsindic,1.,procname) 
    612       else 
    613       CALL Agrif_Interp_variable(q,tabvarsindic) 
    614       CALL Agrif_Bc_variable(q,tabvarsindic,1.) 
     704      CALL Agrif_Interp_variable(q,indic,procname) 
     705      CALL Agrif_Bc_variable(q,indic,1.,procname) 
     706      else 
     707      CALL Agrif_Interp_variable(q,indic) 
     708      CALL Agrif_Bc_variable(q,indic,1.) 
    615709      endif 
    616710 
     
    625719 
    626720      REAL,  DIMENSION(:,:,:,:) :: q 
    627       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     721      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    628722      External :: procname 
    629723      Optional ::  procname 
    630724C 
    631725      if (Agrif_Root()) Return 
     726C       
     727      indic = tabvarsindic 
     728      if (tabvarsindic >=0) then 
     729        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     730          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     731        endif 
     732      endif        
    632733C 
    633734      if (present(procname)) then 
    634       CALL Agrif_Interp_variable(q,tabvarsindic,procname) 
    635       CALL Agrif_Bc_variable(q,tabvarsindic,1.,procname) 
    636       else 
    637       CALL Agrif_Interp_variable(q,tabvarsindic) 
    638       CALL Agrif_Bc_variable(q,tabvarsindic,1.) 
     735      CALL Agrif_Interp_variable(q,indic,procname) 
     736      CALL Agrif_Bc_variable(q,indic,1.,procname) 
     737      else 
     738      CALL Agrif_Interp_variable(q,indic) 
     739      CALL Agrif_Bc_variable(q,indic,1.) 
    639740      endif 
    640741 
     
    798899      External :: procname 
    799900      Optional ::  procname 
    800       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     901      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    801902C         
    802903      REAL, OPTIONAL :: calledweight 
     
    808909C       
    809910      If (Agrif_Root()) Return 
     911C       
     912      indic = tabvarsindic 
     913      if (tabvarsindic >=0) then 
     914        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     915          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     916        endif 
     917      endif          
    810918       
    811919      if ( PRESENT(calledweight) ) then 
     
    817925      endif 
    818926       
    819       if (tabvarsindic <=0) then 
    820       tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) 
     927      if (indic <=0) then 
     928      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
    821929      parenttabvars => tabvars%parent_var 
    822       roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) 
    823       else 
    824       tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) 
    825       parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) 
    826       roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) 
     930      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     931      else 
     932      tabvars=>Agrif_Curgrid % tabvars(indic) 
     933      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     934      roottabvars => Agrif_Mygrid % tabvars(indic) 
    827935      endif 
    828936             
     
    856964      External :: procname 
    857965      Optional ::  procname 
    858       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     966      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    859967C         
    860968      REAL, OPTIONAL :: calledweight 
     
    866974C       
    867975      If (Agrif_Root()) Return 
     976C       
     977      indic = tabvarsindic 
     978      if (tabvarsindic >=0) then 
     979        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     980          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     981        endif 
     982      endif        
    868983       
    869984      if ( PRESENT(calledweight) ) then 
     
    875990      endif 
    876991       
    877       if (tabvarsindic <=0) then 
    878       tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) 
     992      if (indic <=0) then 
     993      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
    879994      parenttabvars => tabvars%parent_var 
    880       roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) 
    881       else 
    882       tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) 
    883       parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) 
    884       roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) 
     995      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     996      else 
     997      tabvars=>Agrif_Curgrid % tabvars(indic) 
     998      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     999      roottabvars => Agrif_Mygrid % tabvars(indic) 
    8851000      endif 
    8861001             
     
    9141029      External :: procname 
    9151030      Optional ::  procname 
    916       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1031      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    9171032C         
    9181033      REAL, OPTIONAL :: calledweight 
     
    9241039C       
    9251040      If (Agrif_Root()) Return 
     1041C       
     1042      indic = tabvarsindic 
     1043      if (tabvarsindic >=0) then 
     1044        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1045          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1046        endif 
     1047      endif        
    9261048       
    9271049      if ( PRESENT(calledweight) ) then 
     
    9331055      endif 
    9341056       
    935       if (tabvarsindic <=0) then 
    936       tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) 
     1057      if (indic <=0) then 
     1058      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
    9371059      parenttabvars => tabvars%parent_var 
    938       roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) 
    939       else 
    940       tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) 
    941       parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) 
    942       roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) 
     1060      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     1061      else 
     1062      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1063      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1064      roottabvars => Agrif_Mygrid % tabvars(indic) 
    9431065      endif 
    9441066             
     
    9721094      External :: procname 
    9731095      Optional ::  procname 
    974       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1096      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    9751097C         
    9761098      REAL, OPTIONAL :: calledweight 
     
    9831105      If (Agrif_Root()) Return 
    9841106       
     1107C       
     1108      indic = tabvarsindic 
     1109      if (tabvarsindic >=0) then 
     1110        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1111          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1112        endif 
     1113      endif       
     1114       
    9851115      if ( PRESENT(calledweight) ) then 
    9861116        weight=calledweight       
     
    9911121      endif 
    9921122       
    993       if (tabvarsindic <=0) then 
    994       tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) 
     1123      if (indic <=0) then 
     1124      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
    9951125      parenttabvars => tabvars%parent_var 
    996       roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) 
    997       else 
    998       tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) 
    999       parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) 
    1000       roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) 
     1126      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     1127      else 
     1128      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1129      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1130      roottabvars => Agrif_Mygrid % tabvars(indic) 
    10011131      endif 
    10021132             
     
    10301160      External :: procname 
    10311161      Optional ::  procname 
    1032       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1162      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    10331163C         
    10341164      REAL, OPTIONAL :: calledweight 
     
    10401170C       
    10411171      If (Agrif_Root()) Return 
     1172C       
     1173      indic = tabvarsindic 
     1174      if (tabvarsindic >=0) then 
     1175        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1176          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1177        endif 
     1178      endif       
    10421179       
    10431180      if ( PRESENT(calledweight) ) then 
     
    10491186      endif 
    10501187       
    1051       if (tabvarsindic <=0) then 
    1052       tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) 
     1188      if (indic <=0) then 
     1189      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
    10531190      parenttabvars => tabvars%parent_var 
    1054       roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) 
    1055       else 
    1056       tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) 
    1057       parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) 
    1058       roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) 
     1191      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     1192      else 
     1193      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1194      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1195      roottabvars => Agrif_Mygrid % tabvars(indic) 
    10591196      endif 
    10601197             
     
    10861223 
    10871224      INTEGER :: tabvarsindic0 ! indice of the variable in tabvars 
    1088       INTEGER :: tabvarsindic  ! indice of the variable in tabvars 
     1225      INTEGER :: tabvarsindic, indic  ! indice of the variable in tabvars 
    10891226      INTEGER :: dimensio  ! indice of the variable in tabvars 
    10901227      External :: procname 
     
    10931230      if (Agrif_Root()) Return 
    10941231C      
    1095       dimensio = Agrif_Mygrid % tabvars(tabvarsindic) % var % nbdim  
     1232      indic = tabvarsindic 
     1233      if (tabvarsindic >=0) then 
     1234        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1235          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1236        endif 
     1237      endif        
     1238C      
     1239      dimensio = Agrif_Mygrid % tabvars(indic) % var % nbdim  
    10961240C 
    10971241      if ( dimensio .EQ. 1 ) then 
    10981242       if (present(procname)) then 
    10991243       Call Agrif_Interp_1D( 
    1100      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1101      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1102      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1244     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1245     & Agrif_Curgrid % parent % tabvars(indic), 
     1246     & Agrif_Curgrid % tabvars(indic), 
    11031247     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array1 ,      
    1104      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1105      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1248     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1249     & Agrif_Mygrid % tabvars(indic) %var % nbdim,procname) 
    11061250       else 
    11071251       Call Agrif_Interp_1D( 
    1108      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1109      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1110      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1252     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1253     & Agrif_Curgrid % parent % tabvars(indic), 
     1254     & Agrif_Curgrid % tabvars(indic), 
    11111255     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array1 ,      
    1112      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1113      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
     1256     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1257     & Agrif_Mygrid % tabvars(indic) %var % nbdim) 
    11141258       endif 
    11151259       endif 
     
    11181262      if (present(procname)) then 
    11191263       Call Agrif_Interp_2D( 
    1120      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1121      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1122      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1264     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1265     & Agrif_Curgrid % parent % tabvars(indic), 
     1266     & Agrif_Curgrid % tabvars(indic), 
    11231267     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array2 ,      
    1124      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1125      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1268     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1269     & Agrif_Mygrid % tabvars(indic) %var % nbdim,procname) 
    11261270      else 
    11271271       Call Agrif_Interp_2D( 
    1128      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1129      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1130      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1272     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1273     & Agrif_Curgrid % parent % tabvars(indic), 
     1274     & Agrif_Curgrid % tabvars(indic), 
    11311275     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array2 ,      
    1132      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1133      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
     1276     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1277     & Agrif_Mygrid % tabvars(indic) %var % nbdim) 
    11341278      endif 
    11351279      endif 
     
    11381282      if (present(procname)) then 
    11391283       Call Agrif_Interp_3D( 
    1140      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1141      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1142      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1284     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1285     & Agrif_Curgrid % parent % tabvars(indic), 
     1286     & Agrif_Curgrid % tabvars(indic), 
    11431287     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array3 ,      
    1144      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1145      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1288     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1289     & Agrif_Mygrid % tabvars(indic) %var % nbdim,procname) 
    11461290      else 
    11471291       Call Agrif_Interp_3D( 
    1148      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1149      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1150      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1292     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1293     & Agrif_Curgrid % parent % tabvars(indic), 
     1294     & Agrif_Curgrid % tabvars(indic), 
    11511295     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array3 ,      
    1152      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1153      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
     1296     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1297     & Agrif_Mygrid % tabvars(indic) %var % nbdim) 
    11541298      endif 
    11551299      endif 
     
    11581302      if (present(procname)) then 
    11591303       Call Agrif_Interp_4D( 
    1160      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1161      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1162      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1304     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1305     & Agrif_Curgrid % parent % tabvars(indic), 
     1306     & Agrif_Curgrid % tabvars(indic), 
    11631307     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array4 ,      
    1164      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1165      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1308     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1309     & Agrif_Mygrid % tabvars(indic) %var % nbdim,procname) 
    11661310      else 
    11671311       Call Agrif_Interp_4D( 
    1168      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1169      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1170      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1312     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1313     & Agrif_Curgrid % parent % tabvars(indic), 
     1314     & Agrif_Curgrid % tabvars(indic), 
    11711315     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array4 ,      
    1172      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1173      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
     1316     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1317     & Agrif_Mygrid % tabvars(indic) %var % nbdim) 
    11741318      endif 
    11751319      endif 
     
    11781322      if (present(procname)) then 
    11791323       Call Agrif_Interp_5D( 
    1180      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1181      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1182      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1324     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1325     & Agrif_Curgrid % parent % tabvars(indic), 
     1326     & Agrif_Curgrid % tabvars(indic), 
    11831327     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array5 ,      
    1184      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1185      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1328     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1329     & Agrif_Mygrid % tabvars(indic) %var % nbdim,procname) 
    11861330      else 
    11871331       Call Agrif_Interp_5D( 
    1188      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1189      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1190      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1332     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1333     & Agrif_Curgrid % parent % tabvars(indic), 
     1334     & Agrif_Curgrid % tabvars(indic), 
    11911335     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array5 ,      
    1192      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1193      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
     1336     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1337     & Agrif_Mygrid % tabvars(indic) %var % nbdim) 
    11941338       endif 
    11951339       endif 
     
    11981342      if (present(procname)) then 
    11991343       Call Agrif_Interp_6D( 
    1200      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1201      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1202      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1344     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1345     & Agrif_Curgrid % parent % tabvars(indic), 
     1346     & Agrif_Curgrid % tabvars(indic), 
    12031347     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array6 ,      
    1204      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1205      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1348     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1349     & Agrif_Mygrid % tabvars(indic) %var % nbdim,procname) 
    12061350      else 
    12071351       Call Agrif_Interp_6D( 
    1208      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1209      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1210      & Agrif_Curgrid % tabvars(tabvarsindic), 
     1352     & Agrif_Mygrid % tabvars(indic) % var %  TypeInterp, 
     1353     & Agrif_Curgrid % parent % tabvars(indic), 
     1354     & Agrif_Curgrid % tabvars(indic), 
    12111355     & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array6 ,      
    1212      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1213      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
     1356     & Agrif_Mygrid % tabvars(indic) % var % restaure, 
     1357     & Agrif_Mygrid % tabvars(indic) %var % nbdim) 
    12141358      endif 
    12151359      endif 
     
    12251369 
    12261370      REAL, DIMENSION(:) :: q 
    1227       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1371      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    12281372      External :: procname 
    12291373      Optional ::  procname 
     1374      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars       
    12301375C 
    12311376      if (Agrif_Root()) Return 
    12321377C       
     1378C       
     1379      indic = tabvarsindic 
     1380      if (tabvarsindic >=0) then 
     1381        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1382          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1383        endif 
     1384      endif  
     1385       
     1386      if (indic <=0) then 
     1387      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     1388      parenttabvars => tabvars%parent_var 
     1389      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     1390      else 
     1391      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1392      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1393      roottabvars => Agrif_Mygrid % tabvars(indic) 
     1394      endif 
     1395       
    12331396      if (present(procname)) then 
    12341397      Call Agrif_Interp_1D( 
    1235      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1236      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1237      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1238      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1239      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1398     & roottabvars % var %  TypeInterp, 
     1399     & parenttabvars, 
     1400     & tabvars,q, 
     1401     & roottabvars % var % restaure, 
     1402     & roottabvars %var % nbdim,procname) 
    12401403      else 
    12411404      Call Agrif_Interp_1D( 
    1242      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1243      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1244      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1245      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1246      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
     1405     & roottabvars % var %  TypeInterp, 
     1406     & parenttabvars, 
     1407     & tabvars,q, 
     1408     & roottabvars % var % restaure, 
     1409     & roottabvars %var % nbdim) 
     1410      
    12471411      endif 
    12481412      Return 
     
    12561420 
    12571421      REAL,  DIMENSION(:,:) :: q 
    1258       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1422      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    12591423      External :: procname 
    12601424      Optional ::  procname 
    1261  
     1425      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars       
    12621426C 
    12631427       if (Agrif_Root()) Return 
    1264 C 
     1428C       
     1429      indic = tabvarsindic 
     1430      if (tabvarsindic >=0) then 
     1431        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1432          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1433        endif 
     1434      endif  
     1435       
     1436      if (indic <=0) then 
     1437      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     1438      parenttabvars => tabvars%parent_var 
     1439      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     1440      if (tabvars%var%restaure) then 
     1441        if (agrif_curgrid%ngridstep == 0) then 
     1442          call AGRIF_CopyFromold_AllOneVar 
     1443     &            (Agrif_Curgrid,Agrif_OldMygrid,indic) 
     1444        endif 
     1445      endif 
     1446      else 
     1447      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1448      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1449      roottabvars => Agrif_Mygrid % tabvars(indic) 
     1450      endif 
     1451 
     1452             
    12651453       if (present(procname)) then 
    12661454       Call Agrif_Interp_2D( 
    1267      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1268      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1269      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1270      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1271      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1455     & roottabvars % var %  TypeInterp, 
     1456     & parenttabvars, 
     1457     & tabvars,q, 
     1458     & roottabvars % var % restaure, 
     1459     & roottabvars %var % nbdim,procname) 
    12721460       else 
    12731461       Call Agrif_Interp_2D( 
    1274      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1275      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1276      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1277      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1278      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
    1279        endif 
     1462     & roottabvars % var %  TypeInterp, 
     1463     & parenttabvars, 
     1464     & tabvars,q, 
     1465     & roottabvars % var % restaure, 
     1466     & roottabvars %var % nbdim) 
     1467      
     1468      endif  
    12801469      Return 
    12811470      End Subroutine Agrif_Interp_var2d 
     
    12881477 
    12891478      REAL,  DIMENSION(:,:,:) :: q 
    1290       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1479      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    12911480      External :: procname 
    12921481      Optional ::  procname 
    1293  
     1482      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars       
    12941483C 
    12951484      if (Agrif_Root()) Return 
    1296 C 
     1485C       
     1486 
     1487      indic = tabvarsindic 
     1488      if (tabvarsindic >=0) then 
     1489        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1490          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1491        endif 
     1492      endif  
     1493       
     1494      if (indic <=0) then 
     1495      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     1496      parenttabvars => tabvars%parent_var 
     1497      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     1498      if (tabvars%var%restaure) then 
     1499        if (agrif_curgrid%ngridstep == 0) then 
     1500          call AGRIF_CopyFromold_AllOneVar 
     1501     &            (Agrif_Curgrid,Agrif_OldMygrid,indic) 
     1502        endif 
     1503      endif       
     1504      else 
     1505      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1506      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1507      roottabvars => Agrif_Mygrid % tabvars(indic) 
     1508      endif 
     1509 
    12971510      if (present(procname)) then 
    12981511      Call Agrif_Interp_3D( 
    1299      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1300      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1301      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1302      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1303      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1512     & roottabvars % var %  TypeInterp, 
     1513     & parenttabvars, 
     1514     & tabvars,q, 
     1515     & roottabvars % var % restaure, 
     1516     & roottabvars %var % nbdim,procname) 
    13041517      else 
    13051518      Call Agrif_Interp_3D( 
    1306      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1307      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1308      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1309      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1310      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
    1311       endif 
     1519     & roottabvars % var %  TypeInterp, 
     1520     & parenttabvars, 
     1521     & tabvars,q, 
     1522     & roottabvars % var % restaure, 
     1523     & roottabvars %var % nbdim) 
     1524      
     1525      endif    
    13121526      Return 
    13131527      End Subroutine Agrif_Interp_var3d 
     
    13201534 
    13211535      REAL,  DIMENSION(:,:,:,:) :: q 
    1322       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1536      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    13231537      External :: procname 
    13241538      Optional ::  procname 
    1325  
     1539      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars   
    13261540C 
    13271541      if (Agrif_Root()) Return 
    1328 C 
     1542C       
     1543      indic = tabvarsindic 
     1544      if (tabvarsindic >=0) then 
     1545        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1546          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1547        endif 
     1548      endif  
     1549       
     1550      if (indic <=0) then 
     1551      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     1552      parenttabvars => tabvars%parent_var 
     1553      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     1554      if (tabvars%var%restaure) then 
     1555        if (agrif_curgrid%ngridstep == 0) then 
     1556          call AGRIF_CopyFromold_AllOneVar 
     1557     &            (Agrif_Curgrid,Agrif_OldMygrid,indic) 
     1558        endif 
     1559      endif       
     1560      else 
     1561      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1562      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1563      roottabvars => Agrif_Mygrid % tabvars(indic) 
     1564      endif 
     1565             
    13291566      if (present(procname)) then 
    13301567      Call Agrif_Interp_4D( 
    1331      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1332      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1333      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1334      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1335      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1568     & roottabvars % var %  TypeInterp, 
     1569     & parenttabvars, 
     1570     & tabvars,q, 
     1571     & roottabvars % var % restaure, 
     1572     & roottabvars %var % nbdim,procname) 
    13361573      else 
    13371574      Call Agrif_Interp_4D( 
    1338      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1339      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1340      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1341      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1342      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
    1343       endif 
     1575     & roottabvars % var %  TypeInterp, 
     1576     & parenttabvars, 
     1577     & tabvars,q, 
     1578     & roottabvars % var % restaure, 
     1579     & roottabvars %var % nbdim) 
     1580      
     1581      endif       
     1582 
    13441583      Return 
    13451584      End Subroutine Agrif_Interp_var4d      
     
    13521591 
    13531592      REAL,  DIMENSION(:,:,:,:,:) :: q 
    1354       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1593      INTEGER :: tabvarsindic, indic ! indice of the variable in tabvars 
    13551594      External :: procname 
    13561595      Optional ::  procname 
    1357  
     1596      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars       
    13581597C 
    13591598      if (Agrif_Root()) Return 
    1360 C 
     1599C       
     1600 
     1601      indic = tabvarsindic 
     1602      if (tabvarsindic >=0) then 
     1603        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1604          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1605        endif 
     1606      endif  
     1607       
     1608      if (indic <=0) then 
     1609      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     1610      parenttabvars => tabvars%parent_var 
     1611      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     1612      else 
     1613      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1614      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1615      roottabvars => Agrif_Mygrid % tabvars(indic) 
     1616      endif 
     1617       
    13611618      if (present(procname)) then 
    13621619      Call Agrif_Interp_5D( 
    1363      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1364      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1365      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1366      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1367      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) 
     1620     & roottabvars % var %  TypeInterp, 
     1621     & parenttabvars, 
     1622     & tabvars,q, 
     1623     & roottabvars % var % restaure, 
     1624     & roottabvars %var % nbdim,procname) 
    13681625      else 
    13691626      Call Agrif_Interp_5D( 
    1370      & Agrif_Mygrid % tabvars(tabvarsindic) % var %  TypeInterp, 
    1371      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1372      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1373      & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, 
    1374      & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) 
     1627     & roottabvars % var %  TypeInterp, 
     1628     & parenttabvars, 
     1629     & tabvars,q, 
     1630     & roottabvars % var % restaure, 
     1631     & roottabvars %var % nbdim) 
     1632      
    13751633      endif 
    13761634      Return 
     
    15381796 
    15391797      REAL,  DIMENSION(:) :: q 
    1540       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1798      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
    15411799      External :: procname 
    15421800      Optional ::  procname       
     
    15441802      INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 
    15451803      INTEGER, DIMENSION(2), OPTIONAL :: locupdate2         
     1804      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars    
    15461805C       
    15471806      if (Agrif_Root()) Return 
    15481807C      
     1808 
     1809      indic = tabvarsindic 
     1810      if (tabvarsindic >=0) then 
     1811        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1812          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1813        endif 
     1814      endif  
     1815       
     1816      if (indic <=0) then 
     1817      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     1818      parenttabvars => tabvars%parent_var 
     1819      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic)       
     1820      else 
     1821      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1822      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1823      roottabvars => Agrif_Mygrid % tabvars(indic) 
     1824      endif 
     1825       
    15491826      IF (present(locupdate)) THEN 
    1550       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:1)  
     1827      tabvars%var % updateinf(1:1)  
    15511828     &      = locupdate(1) 
    1552       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:1)  
     1829      tabvars%var % updatesup(1:1)  
    15531830     &      = locupdate(2) 
    15541831      ELSE 
    1555       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:1)  
     1832      tabvars%var % updateinf(1:1)  
    15561833     &      = -99 
    1557       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:1)  
     1834      tabvars%var % updatesup(1:1)  
    15581835     &      = -99 
    15591836      ENDIF 
    15601837       
    15611838      IF (present(locupdate1)) THEN 
    1562       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1)  
     1839      tabvars%var % updateinf(1)  
    15631840     &      = locupdate1(1) 
    1564       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1)  
     1841      tabvars%var % updatesup(1)  
    15651842     &      = locupdate1(2) 
    15661843      ENDIF   
    15671844       
    15681845      IF (present(locupdate2)) THEN 
    1569       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2)  
     1846      tabvars%var % updateinf(2)  
    15701847     &      = locupdate2(1) 
    1571       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2)  
     1848      tabvars%var % updatesup(2)  
    15721849     &      = locupdate2(2) 
    15731850      ENDIF        
     
    15751852      IF (present(procname)) THEN 
    15761853      Call Agrif_Update_1D( 
    1577      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1578      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1579      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1580      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1581      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, 
     1854     & roottabvars % var % typeupdate, 
     1855     & parenttabvars, 
     1856     & tabvars,q, 
     1857     & tabvars % var % updateinf, 
     1858     & tabvars % var % updatesup, 
    15821859     & procname) 
    1583       ELSE 
     1860      ELSE  
    15841861      Call Agrif_Update_1D( 
    1585      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1586      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1587      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1588      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1589      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup)       
     1862     & roottabvars % var % typeupdate, 
     1863     & parenttabvars, 
     1864     & tabvars,q, 
     1865     & tabvars % var % updateinf, 
     1866     & tabvars % var % updatesup)        
    15901867      ENDIF 
    15911868 
     
    16071884      INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 
    16081885      INTEGER, DIMENSION(2), OPTIONAL :: locupdate2         
    1609       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1886      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
     1887      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars        
    16101888C       
    16111889      IF (Agrif_Root()) RETURN 
    16121890       
    16131891C  
     1892      indic = tabvarsindic 
     1893      if (tabvarsindic >=0) then 
     1894        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1895          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1896        endif 
     1897      endif  
     1898       
     1899      if (indic <=0) then 
     1900      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     1901      parenttabvars => tabvars%parent_var 
     1902      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic)       
     1903      else 
     1904      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1905      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1906      roottabvars => Agrif_Mygrid % tabvars(indic) 
     1907      endif 
     1908       
    16141909      IF (present(locupdate)) THEN 
    1615       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:2)  
     1910      tabvars%var % updateinf(1:2)  
    16161911     &      = locupdate(1) 
    1617       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:2)  
     1912      tabvars%var % updatesup(1:2)  
    16181913     &      = locupdate(2) 
    16191914      ELSE 
    1620       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:2)  
     1915      tabvars%var % updateinf(1:2)  
    16211916     &      = -99 
    1622       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:2)  
     1917      tabvars%var % updatesup(1:2)  
    16231918     &      = -99 
    16241919      ENDIF 
    16251920       
    16261921      IF (present(locupdate1)) THEN 
    1627       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1)  
     1922      tabvars%var % updateinf(1)  
    16281923     &      = locupdate1(1) 
    1629       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1)  
     1924      tabvars%var % updatesup(1)  
    16301925     &      = locupdate1(2) 
    16311926      ENDIF   
    16321927       
    16331928      IF (present(locupdate2)) THEN 
    1634       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2)  
     1929      tabvars%var % updateinf(2)  
    16351930     &      = locupdate2(1) 
    1636       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2)  
     1931      tabvars%var % updatesup(2)  
    16371932     &      = locupdate2(2) 
    16381933      ENDIF 
     
    16401935      IF (present(procname)) THEN 
    16411936      Call Agrif_Update_2D( 
    1642      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1643      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1644      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1645      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1646      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, 
     1937     & roottabvars % var % typeupdate, 
     1938     & parenttabvars, 
     1939     & tabvars,q, 
     1940     & tabvars % var % updateinf, 
     1941     & tabvars % var % updatesup, 
    16471942     & procname) 
    1648       ELSE 
     1943      ELSE  
    16491944      Call Agrif_Update_2D( 
    1650      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1651      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1652      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1653      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1654      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup)       
     1945     & roottabvars % var % typeupdate, 
     1946     & parenttabvars, 
     1947     & tabvars,q, 
     1948     & tabvars % var % updateinf, 
     1949     & tabvars % var % updatesup)        
    16551950      ENDIF 
    16561951 
     
    16721967      INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 
    16731968      INTEGER, DIMENSION(2), OPTIONAL :: locupdate2         
    1674       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     1969      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
     1970      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars        
    16751971C       
    16761972      IF (Agrif_Root()) RETURN  
    16771973C       
    1678  
     1974      indic = tabvarsindic 
     1975      if (tabvarsindic >=0) then 
     1976        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     1977          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     1978        endif 
     1979      endif  
     1980       
     1981      if (indic <=0) then 
     1982      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     1983      parenttabvars => tabvars%parent_var 
     1984      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic)       
     1985      else 
     1986      tabvars=>Agrif_Curgrid % tabvars(indic) 
     1987      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     1988      roottabvars => Agrif_Mygrid % tabvars(indic) 
     1989      endif 
     1990       
    16791991      IF (present(locupdate)) THEN 
    1680       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:3)  
     1992      tabvars%var % updateinf(1:3)  
    16811993     &      = locupdate(1) 
    1682       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:3)  
     1994      tabvars%var % updatesup(1:3)  
    16831995     &      = locupdate(2) 
    16841996      ELSE 
    1685       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:3)  
     1997      tabvars%var % updateinf(1:3)  
    16861998     &      = -99 
    1687       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:3)  
     1999      tabvars%var % updatesup(1:3)  
    16882000     &      = -99 
    16892001      ENDIF       
    16902002       
    16912003      IF (present(locupdate1)) THEN 
    1692       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1)  
     2004      tabvars%var % updateinf(1)  
    16932005     &      = locupdate1(1) 
    1694       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1)  
     2006      tabvars%var % updatesup(1)  
    16952007     &      = locupdate1(2) 
    16962008      ENDIF   
    16972009       
    16982010      IF (present(locupdate2)) THEN 
    1699       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2)  
     2011      tabvars%var % updateinf(2)  
    17002012     &      = locupdate2(1) 
    1701       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2)  
     2013      tabvars%var % updatesup(2)  
    17022014     &      = locupdate2(2) 
    17032015      ENDIF 
     
    17052017      IF (present(procname)) THEN 
    17062018      Call Agrif_Update_3D( 
    1707      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1708      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1709      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1710      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1711      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, 
     2019     & roottabvars % var % typeupdate, 
     2020     & parenttabvars, 
     2021     & tabvars,q, 
     2022     & tabvars % var % updateinf, 
     2023     & tabvars % var % updatesup, 
    17122024     & procname) 
    1713       ELSE 
     2025      ELSE  
    17142026      Call Agrif_Update_3D( 
    1715      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1716      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1717      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1718      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1719      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup)       
     2027     & roottabvars % var % typeupdate, 
     2028     & parenttabvars, 
     2029     & tabvars,q, 
     2030     & tabvars % var % updateinf, 
     2031     & tabvars % var % updatesup)        
    17202032      ENDIF 
    17212033 
     
    17372049      INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 
    17382050      INTEGER, DIMENSION(2), OPTIONAL :: locupdate2         
    1739       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     2051      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
     2052      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars        
    17402053C       
    17412054      IF (Agrif_Root()) RETURN 
     2055      indic = tabvarsindic 
     2056      if (tabvarsindic >=0) then 
     2057        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     2058          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     2059        endif 
     2060      endif  
     2061       
     2062      if (indic <=0) then 
     2063      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     2064      parenttabvars => tabvars%parent_var 
     2065      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic)      
     2066      else 
     2067      tabvars=>Agrif_Curgrid % tabvars(indic) 
     2068      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     2069      roottabvars => Agrif_Mygrid % tabvars(indic) 
     2070      endif       
    17422071C       
    17432072      IF (present(locupdate)) THEN 
    1744       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:4)  
     2073      tabvars%var % updateinf(1:4)  
    17452074     &      = locupdate(1) 
    1746       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:4)  
     2075      tabvars%var % updatesup(1:4)  
    17472076     &      = locupdate(2) 
    17482077      ELSE 
    1749       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:4)  
     2078      tabvars%var % updateinf(1:4)  
    17502079     &      = -99 
    1751       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:4)  
     2080      tabvars%var % updatesup(1:4)  
    17522081     &      = -99 
    17532082      ENDIF 
    17542083       
    17552084      IF (present(locupdate1)) THEN 
    1756       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1)  
     2085      tabvars%var % updateinf(1)  
    17572086     &      = locupdate1(1) 
    1758       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1)  
     2087      tabvars%var % updatesup(1)  
    17592088     &      = locupdate1(2) 
    17602089      ENDIF   
    17612090       
    17622091      IF (present(locupdate2)) THEN 
    1763       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2)  
     2092      tabvars%var % updateinf(2)  
    17642093     &      = locupdate2(1) 
    1765       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2)  
     2094      tabvars%var % updatesup(2)  
    17662095     &      = locupdate2(2) 
    17672096      ENDIF 
     
    17692098      IF (present(procname)) THEN 
    17702099      Call Agrif_Update_4D( 
    1771      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1772      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1773      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1774      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1775      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, 
     2100     & roottabvars % var % typeupdate, 
     2101     & parenttabvars, 
     2102     & tabvars,q, 
     2103     & tabvars % var % updateinf, 
     2104     & tabvars % var % updatesup, 
    17762105     & procname) 
    1777       ELSE 
     2106      ELSE  
    17782107      Call Agrif_Update_4D( 
    1779      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1780      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1781      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1782      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1783      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup)       
     2108     & roottabvars % var % typeupdate, 
     2109     & parenttabvars, 
     2110     & tabvars,q, 
     2111     & tabvars % var % updateinf, 
     2112     & tabvars % var % updatesup)        
    17842113      ENDIF 
    17852114 
     
    18012130      INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 
    18022131      INTEGER, DIMENSION(2), OPTIONAL :: locupdate2         
    1803       INTEGER :: tabvarsindic ! indice of the variable in tabvars 
     2132      INTEGER :: tabvarsindic,indic ! indice of the variable in tabvars 
     2133      TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars        
    18042134C 
    18052135      IF (Agrif_Root()) RETURN 
    18062136C       
     2137      indic = tabvarsindic 
     2138      if (tabvarsindic >=0) then 
     2139        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     2140          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     2141        endif 
     2142      endif  
     2143       
     2144      if (indic <=0) then 
     2145      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     2146      parenttabvars => tabvars%parent_var 
     2147      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic)       
     2148      else 
     2149      tabvars=>Agrif_Curgrid % tabvars(indic) 
     2150      parenttabvars => Agrif_Curgrid % parent % tabvars(indic) 
     2151      roottabvars => Agrif_Mygrid % tabvars(indic) 
     2152      endif 
     2153        
    18072154      IF (present(locupdate)) THEN 
    1808       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:5)  
     2155      tabvars%var % updateinf(1:5)  
    18092156     &      = locupdate(1) 
    1810       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:5)  
     2157      tabvars%var % updatesup(1:5)  
    18112158     &      = locupdate(2) 
    18122159      ELSE 
    1813       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:5)  
     2160      tabvars%var % updateinf(1:5)  
    18142161     &      = -99 
    1815       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:5)  
     2162      tabvars%var % updatesup(1:5)  
    18162163     &      = -99 
    18172164      ENDIF 
    18182165       
    18192166      IF (present(locupdate1)) THEN 
    1820       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1)  
     2167      tabvars%var % updateinf(1)  
    18212168     &      = locupdate1(1) 
    1822       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1)  
     2169      tabvars%var % updatesup(1)  
    18232170     &      = locupdate1(2) 
    18242171      ENDIF   
    18252172       
    18262173      IF (present(locupdate2)) THEN 
    1827       Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2)  
     2174      tabvars%var % updateinf(2)  
    18282175     &      = locupdate2(1) 
    1829       Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2)  
     2176      tabvars%var % updatesup(2)  
    18302177     &      = locupdate2(2) 
    18312178      ENDIF 
     
    18332180      IF (present(procname)) THEN 
    18342181      Call Agrif_Update_5D( 
    1835      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1836      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1837      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1838      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1839      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, 
     2182     & roottabvars % var % typeupdate, 
     2183     & parenttabvars, 
     2184     & tabvars,q, 
     2185     & tabvars % var % updateinf, 
     2186     & tabvars % var % updatesup, 
    18402187     & procname) 
    1841       ELSE 
     2188      ELSE  
    18422189      Call Agrif_Update_5D( 
    1843      & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, 
    1844      & Agrif_Curgrid % parent % tabvars(tabvarsindic), 
    1845      & Agrif_Curgrid % tabvars(tabvarsindic),q, 
    1846      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, 
    1847      & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup)       
     2190     & roottabvars % var % typeupdate, 
     2191     & parenttabvars, 
     2192     & tabvars,q, 
     2193     & tabvars % var % updateinf, 
     2194     & tabvars % var % updatesup)        
    18482195      ENDIF 
    18492196 
     
    19522299      End Subroutine Agrif_Flux_Correction 
    19532300 
    1954       Subroutine Agrif_Declare_Variable(posvar,firstpoint, 
    1955      &    raf,lb,ub,varid) 
    1956       character*(80) :: variablename 
    1957       Type(Agrif_List_Variables), Pointer :: newvariable,newvariablep 
    1958       INTEGER, DIMENSION(:) :: posvar 
    1959       INTEGER, DIMENSION(:) :: lb,ub 
    1960       INTEGER, DIMENSION(:) :: firstpoint 
    1961       CHARACTER(*) ,DIMENSION(:) :: raf         
    1962       TYPE(Agrif_Pvariable), Pointer :: parent_var,root_var 
    1963       INTEGER :: dimensio 
    1964       INTEGER :: varid 
    1965              
    1966       if (agrif_root()) return 
    1967  
    1968       dimensio = SIZE(posvar) 
    1969 C 
    1970 C     
    1971       Allocate(newvariable) 
    1972       Allocate(newvariable%pvar) 
    1973       Allocate(newvariable%pvar%var) 
    1974       Allocate(newvariable%pvar%var%posvar(dimensio)) 
    1975       Allocate(newvariable%pvar%var%interptab(dimensio)) 
    1976       newvariable%pvar%var%variablename = variablename 
    1977       newvariable%pvar%var%interptab = raf 
    1978       newvariable%pvar%var%nbdim = dimensio 
    1979       newvariable%pvar%var%posvar = posvar 
    1980       newvariable%pvar%var%point(1:dimensio) = firstpoint 
    1981       newvariable%pvar%var%lb(1:dimensio) = lb(1:dimensio) 
    1982       newvariable%pvar%var%ub(1:dimensio) = ub(1:dimensio) 
    1983        
    1984       newvariable % nextvariable => Agrif_Curgrid%variables 
    1985        
    1986       Agrif_Curgrid%variables => newvariable 
    1987       Agrif_Curgrid%Nbvariables = Agrif_Curgrid%Nbvariables + 1 
    1988        
    1989       varid = -Agrif_Curgrid%Nbvariables 
    1990        
    1991        if (agrif_curgrid%parent%nbvariables < agrif_curgrid%nbvariables) 
    1992      &       then 
    1993       Allocate(newvariablep) 
    1994       Allocate(newvariablep%pvar) 
    1995       Allocate(newvariablep%pvar%var)       
    1996       Allocate(newvariablep%pvar%var%posvar(dimensio)) 
    1997       Allocate(newvariablep%pvar%var%interptab(dimensio)) 
    1998       newvariablep%pvar%var%variablename = variablename 
    1999       newvariablep%pvar%var%interptab = raf 
    2000       newvariablep%pvar%var%nbdim = dimensio 
    2001       newvariablep%pvar%var%posvar = posvar 
    2002       newvariablep%pvar%var%point(1:dimensio) = firstpoint 
    2003        
    2004       newvariablep % nextvariable => Agrif_Curgrid%parent%variables 
    2005        
    2006       Agrif_Curgrid%parent%variables => newvariablep        
    2007        
    2008       Agrif_Curgrid%parent%Nbvariables =  
    2009      &    Agrif_Curgrid%parent%Nbvariables + 1 
    2010       parent_var=>newvariablep%pvar 
    2011       else 
    2012       parent_var=>Agrif_Search_Variable 
    2013      &              (Agrif_Curgrid%parent,Agrif_Curgrid%nbvariables) 
    2014        endif 
    2015         
    2016        newvariable%pvar%parent_var=>parent_var 
    2017        
    2018       root_var=>Agrif_Search_Variable 
    2019      &              (Agrif_Mygrid,Agrif_Curgrid%nbvariables) 
    2020       
    2021       newvariable%pvar%var%root_var=>root_var%var 
    2022        
    2023             
    2024       End Subroutine Agrif_Declare_Variable 
    2025  
    2026       FUNCTION Agrif_Search_Variable(grid,varid) 
    2027       integer :: varid 
    2028       Type(Agrif_Pvariable), Pointer :: Agrif_Search_variable 
    2029       Type(Agrif_grid), Pointer :: grid 
    2030        
    2031       Type(Agrif_List_Variables), pointer :: parcours 
    2032       Logical :: foundvariable 
    2033       integer nb 
    2034        
    2035       foundvariable = .FALSE. 
    2036       parcours => grid%variables 
    2037        
    2038       do nb=1,varid-1 
    2039          parcours => parcours%nextvariable 
    2040       End Do 
    2041        
    2042       Agrif_Search_variable => parcours%pvar 
    2043        
    2044        
    2045       End Function Agrif_Search_variable 
    2046                                
     2301 
     2302       
    20472303      Subroutine Agrif_Declare_Profile_flux(profilename,posvar, 
    20482304     &    firstpoint,raf) 
     
    20712327       
    20722328      End Subroutine Agrif_Declare_Profile_flux 
     2329 
     2330      Subroutine Agrif_Save_ForRestore0D(tabvarsindic0,tabvarsindic) 
     2331      integer :: tabvarsindic0, tabvarsindic 
     2332      integer :: dimensio 
    20732333               
     2334      dimensio =  Agrif_Mygrid % tabvars(tabvarsindic0) % var % nbdim 
     2335       
     2336      select case(dimensio) 
     2337      case(2) 
     2338          call Agrif_Save_ForRestore2D( 
     2339     &      Agrif_Curgrid % tabvars(tabvarsindic0) % var % array2, 
     2340     &      tabvarsindic) 
     2341      case(3) 
     2342          call Agrif_Save_ForRestore3D( 
     2343     &      Agrif_Curgrid % tabvars(tabvarsindic0) % var % array3, 
     2344     &      tabvarsindic) 
     2345      case(4) 
     2346          call Agrif_Save_ForRestore4D( 
     2347     &      Agrif_Curgrid % tabvars(tabvarsindic0) % var % array4, 
     2348     &      tabvarsindic)      
     2349      end select 
     2350 
     2351      Return 
     2352      End Subroutine Agrif_Save_ForRestore0D  
     2353       
     2354     
     2355       
     2356      Subroutine Agrif_Save_ForRestore2D(q,tabvarsindic) 
     2357      real,dimension(:,:) :: q 
     2358      integer :: tabvarsindic, indic 
     2359      TYPE(Agrif_PVariable),Pointer ::tabvars, roottabvars     
     2360  
     2361       indic = tabvarsindic 
     2362      if (tabvarsindic >=0) then 
     2363        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     2364          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     2365        endif 
     2366      endif 
     2367       
     2368      if (indic <=0) then 
     2369      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     2370      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     2371      else 
     2372      tabvars=>Agrif_Curgrid % tabvars(indic) 
     2373      roottabvars => Agrif_Mygrid % tabvars(indic) 
     2374      endif       
     2375      if (.not.allocated(tabvars%var%array2)) then       
     2376      allocate(tabvars%var%array2(tabvars%var%lb(1):tabvars%var%ub(1), 
     2377     &                            tabvars%var%lb(2):tabvars%var%ub(2))) 
     2378      endif 
     2379      tabvars%var%array2 = q 
     2380      roottabvars%var%restaure = .true. 
     2381       
     2382      Return 
     2383      End Subroutine Agrif_Save_ForRestore2D   
     2384       
     2385      Subroutine Agrif_Save_ForRestore3D(q,tabvarsindic) 
     2386      real,dimension(:,:,:) :: q 
     2387      integer :: tabvarsindic, indic 
     2388      TYPE(Agrif_PVariable),Pointer ::tabvars, roottabvars     
     2389  
     2390       indic = tabvarsindic 
     2391      if (tabvarsindic >=0) then 
     2392        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     2393          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     2394        endif 
     2395      endif 
     2396       
     2397      if (indic <=0) then 
     2398      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     2399      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     2400      else 
     2401      tabvars=>Agrif_Curgrid % tabvars(indic) 
     2402      roottabvars => Agrif_Mygrid % tabvars(indic) 
     2403      endif       
     2404 
     2405      if (.not.allocated(tabvars%var%array3)) then 
     2406      allocate(tabvars%var%array3(tabvars%var%lb(1):tabvars%var%ub(1), 
     2407     &                            tabvars%var%lb(2):tabvars%var%ub(2), 
     2408     &                            tabvars%var%lb(3):tabvars%var%ub(3))) 
     2409      endif 
     2410      tabvars%var%array3 = q 
     2411      roottabvars%var%restaure = .true. 
     2412       
     2413      Return 
     2414      End Subroutine Agrif_Save_ForRestore3D 
     2415       
     2416      Subroutine Agrif_Save_ForRestore4D(q,tabvarsindic) 
     2417      real,dimension(:,:,:,:) :: q 
     2418      integer :: tabvarsindic, indic 
     2419      TYPE(Agrif_PVariable),Pointer ::tabvars, roottabvars     
     2420  
     2421       indic = tabvarsindic 
     2422      if (tabvarsindic >=0) then 
     2423        if (agrif_curgrid%tabvars(tabvarsindic)%var%nbdim == 0) then 
     2424          indic = agrif_curgrid%tabvars(tabvarsindic)%var%iarray0 
     2425        endif 
     2426      endif 
     2427       
     2428      if (indic <=0) then 
     2429      tabvars => Agrif_Search_Variable(Agrif_Curgrid,-indic) 
     2430      roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-indic) 
     2431      else 
     2432      tabvars=>Agrif_Curgrid % tabvars(indic) 
     2433      roottabvars => Agrif_Mygrid % tabvars(indic) 
     2434      endif 
     2435 
     2436      if (.not.allocated(tabvars%var%array4)) then          
     2437      allocate(tabvars%var%array4(tabvars%var%lb(1):tabvars%var%ub(1), 
     2438     &                            tabvars%var%lb(2):tabvars%var%ub(2), 
     2439     &                            tabvars%var%lb(3):tabvars%var%ub(3), 
     2440     &                            tabvars%var%lb(4):tabvars%var%ub(4))) 
     2441      endif 
     2442      tabvars%var%array4 = q 
     2443      roottabvars%var%restaure = .true. 
     2444       
     2445      Return 
     2446      End Subroutine Agrif_Save_ForRestore4D                    
    20742447C 
    20752448      End module Agrif_bcfunction 
Note: See TracChangeset for help on using the changeset viewer.