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.
modsauv.F in branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modsauv.F @ 10251

Last change on this file since 10251 was 10251, checked in by kingr, 5 years ago

Rolled back to r10247 - i.e., undid merge of pkg br and 3.6_stable br

File size: 26.2 KB
Line 
1!
2! $Id: modsauv.F 2715 2011-03-30 15:58:35Z rblod $
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_Save 
26C
27      Module Agrif_Save 
28C
29CCC   Description:   
30CCC   Module for the initialization by copy of the grids created by clustering.
31C
32C     Modules used:
33C     
34      Use Agrif_Types
35      Use Agrif_Link
36      Use Agrif_Arrays
37      Use Agrif_Variables
38C
39      IMPLICIT NONE
40C
41      Contains
42C     Definition of procedures contained in this module.
43C
44C
45      Subroutine Agrif_Deallocate_Arrays(Var)
46      type(Agrif_Variable), pointer :: Var
47     
48            if (ALLOCATED(var%array1)) then
49               Deallocate(var%array1)
50            endif
51            if (ALLOCATED(var%array2)) then
52               Deallocate(var%array2)
53            endif
54            if (ALLOCATED(var%array3)) then
55               Deallocate(var%array3)
56            endif
57            if (ALLOCATED(var%array4)) then
58               Deallocate(var%array4)
59            endif
60            if (ALLOCATED(var%array5)) then
61               Deallocate(var%array5)
62            endif
63            if (ALLOCATED(var%array6)) then
64               Deallocate(var%array6)
65            endif
66C
67            if (ALLOCATED(var%darray1)) then
68               Deallocate(var%darray1)
69            endif
70            if (ALLOCATED(var%darray2)) then
71               Deallocate(var%darray2)
72            endif
73            if (ALLOCATED(var%darray3)) then
74               Deallocate(var%darray3)
75            endif
76            if (ALLOCATED(var%darray4)) then
77               Deallocate(var%darray4)
78            endif
79            if (ALLOCATED(var%darray5)) then
80               Deallocate(var%darray5)
81            endif
82            if (ALLOCATED(var%darray6)) then
83               Deallocate(var%darray6)
84            endif
85C
86            if (ALLOCATED(var%larray1)) then
87               Deallocate(var%larray1)
88            endif
89            if (ALLOCATED(var%larray2)) then
90               Deallocate(var%larray2)
91            endif
92            if (ALLOCATED(var%larray3)) then
93               Deallocate(var%larray3)
94            endif
95            if (ALLOCATED(var%larray4)) then
96               Deallocate(var%larray4)
97            endif
98            if (ALLOCATED(var%larray5)) then
99               Deallocate(var%larray5)
100            endif
101            if (ALLOCATED(var%larray6)) then
102               Deallocate(var%larray6)
103            endif
104C
105            if (ALLOCATED(var%iarray1)) then
106               Deallocate(var%iarray1)
107            endif
108            if (ALLOCATED(var%iarray2)) then
109               Deallocate(var%iarray2)
110            endif
111            if (ALLOCATED(var%iarray3)) then
112               Deallocate(var%iarray3)
113            endif
114            if (ALLOCATED(var%iarray4)) then
115               Deallocate(var%iarray4)
116            endif
117            if (ALLOCATED(var%iarray5)) then
118               Deallocate(var%iarray5)
119            endif
120            if (ALLOCATED(var%iarray6)) then
121               Deallocate(var%iarray6)
122            endif
123C
124            if (ALLOCATED(var%carray1)) then
125               Deallocate(var%carray1)
126            endif
127            if (ALLOCATED(var%carray2)) then
128               Deallocate(var%carray2)
129            endif
130C
131            if (associated(var%oldvalues2D)) then
132               Deallocate(var%oldvalues2D)
133            endif
134            if (associated(var%interpIndex)) then
135               Deallocate(var%interpIndex)
136            endif
137           
138            if (associated(var%posvar)) then
139               Deallocate(var%posvar)
140            endif     
141           
142            if (associated(var%interptab)) then
143               Deallocate(var%interptab)
144            endif
145           
146      Return
147      End Subroutine Agrif_Deallocate_Arrays   
148C
149C     **************************************************************************
150CCC   Subroutine Agrif_Free_data_before
151C     **************************************************************************
152C
153      Subroutine Agrif_Free_data_before(Agrif_Gr)
154C
155CCC   Description:
156CCC   
157C
158CC    Method:
159CC           
160C
161C     Declarations:
162C
163     
164C     
165C     Pointer argument   
166      TYPE(Agrif_Grid),pointer   :: Agrif_Gr ! Pointer on the current grid
167      INTEGER i 
168      Type(Agrif_List_Variables), pointer :: parcours     
169C
170C   
171      do i = 1 , AGRIF_NbVariables
172         if ( .NOT. Agrif_Mygrid % tabvars(i) % var % restaure) then
173C 
174            call Agrif_Deallocate_Arrays(Agrif_Gr%tabvars(i)%var)
175           
176       endif
177           
178C
179       if (associated(Agrif_Gr%tabvars(i)%var%list_interp)) then
180         Call Agrif_Free_list_interp
181     &                          (Agrif_Gr%tabvars(i)%var%list_interp)
182       endif                             
183C
184       if ( .NOT. Agrif_Mygrid % tabvars(i) % var % restaure) then
185            Deallocate(Agrif_Gr%tabvars(i)%var)
186C 
187        endif
188      enddo
189     
190      parcours => Agrif_Gr%variables
191   
192      do i=1,Agrif_Gr%NbVariables
193       if (.NOT. parcours%pvar%var%root_var%restaure) then
194        call Agrif_Deallocate_Arrays(parcours%pvar%var)
195       endif
196        if (associated(parcours%pvar%var%list_interp)) then
197           Call Agrif_Free_list_interp
198     &                          (parcours%pvar%var%list_interp)
199         endif                             
200C
201       if ( .NOT. parcours%pvar%var%root_var % restaure) then
202            Deallocate(parcours%pvar%var)
203C 
204        endif
205        parcours => parcours%nextvariable
206      enddo
207C
208C 
209C
210      if ( Agrif_USE_ONLY_FIXED_GRIDS .EQ. 0 ) then
211         if ( Agrif_Probdim .EQ. 1 ) Deallocate(Agrif_Gr%tabpoint1D)
212         if ( Agrif_Probdim .EQ. 2 ) Deallocate(Agrif_Gr%tabpoint2D)
213         if ( Agrif_Probdim .EQ. 3 ) Deallocate(Agrif_Gr%tabpoint3D)
214      endif
215C         
216      Return
217C
218C
219      End Subroutine Agrif_Free_data_before
220C
221C
222      Recursive Subroutine Agrif_Free_list_interp(list_interp)
223      TYPE(Agrif_List_Interp_Loc), Pointer :: list_interp
224     
225      if (associated(list_interp%suiv))
226     &      Call Agrif_Free_list_interp(list_interp%suiv)
227     
228#ifdef key_mpp_mpi
229       Deallocate(list_interp%interp_loc%tab4t)
230       Deallocate(list_interp%interp_loc%memberinall)
231       Deallocate(list_interp%interp_loc%sendtoproc1)
232       Deallocate(list_interp%interp_loc%recvfromproc1)
233#endif
234       Deallocate(list_interp%interp_loc)
235       Deallocate(list_interp)
236       Nullify(list_interp)
237
238      End Subroutine Agrif_Free_list_interp     
239C
240C     **************************************************************************
241CCC   Subroutine Agrif_Free_data_after
242C     **************************************************************************
243C     
244      Subroutine Agrif_Free_data_after(Agrif_Gr)
245C
246CCC   Description:
247CCC   
248C
249CC    Method:
250CC           
251C
252C     Declarations:
253C
254C
255     
256C
257C     Pointer argument   
258      TYPE(Agrif_Grid),pointer   :: Agrif_Gr  ! Pointer on the current grid
259      INTEGER i
260      Type(Agrif_List_Variables), pointer :: parcours, rootparcours
261C
262C     
263      do i = 1 , AGRIF_NbVariables
264         if ( Agrif_Mygrid % tabvars(i) % var % restaure) then
265           
266         call Agrif_Deallocate_Arrays(Agrif_Gr%tabvars(i)%var)
267         !
268            Deallocate(Agrif_Gr%tabvars(i)%var)
269!
270         endif
271      enddo
272     
273      parcours => Agrif_Gr%variables
274      rootparcours=>Agrif_Mygrid%variables
275   
276      do i=1,Agrif_Gr%NbVariables
277       if (rootparcours%pvar%var%restaure) then
278        call Agrif_Deallocate_Arrays(parcours%pvar%var)
279     
280            Deallocate(parcours%pvar%var)
281C 
282        endif
283        parcours => parcours%nextvariable
284        rootparcours => rootparcours%nextvariable       
285      enddo
286           
287C
288C     
289      Return
290C
291C
292      End Subroutine Agrif_Free_data_after     
293C
294C
295CC    **************************************************************************
296CCC   Subroutine AGRIF_CopyFromold_All
297C     **************************************************************************
298C
299      Recursive Subroutine AGRIF_CopyFromold_All(g,oldchildgrids)
300C
301CCC   Description:
302CCC   Routine called in the Agrif_Init_Hierarchy procedure 
303C       (Agrif_Clustering module). 
304C
305CC    Method:       
306C
307C     Declarations:
308C
309     
310C     
311C     Pointer argument   
312      TYPE(AGRIF_grid),pointer   :: g ! Pointer on the current grid
313      TYPE(AGRIF_pgrid),pointer   :: oldchildgrids
314C
315C     Local pointer
316      TYPE(AGRIF_pgrid),pointer  :: parcours ! Pointer for the recursive
317                                             ! procedure
318      REAL g_eps,eps,oldgrid_eps
319      INTEGER :: out
320      INTEGER :: iii
321C
322      out = 0
323C                                                                           
324      parcours => oldchildgrids 
325C
326      do while (associated(parcours))
327C 
328        if ((.NOT. g % fixed) .AND. (parcours % gr %oldgrid)) then
329C       
330            g_eps = huge(1.)
331            oldgrid_eps = huge(1.)
332            do iii = 1 , Agrif_Probdim
333               g_eps = min(g_eps,g % Agrif_d(iii))
334               oldgrid_eps = min(oldgrid_eps,
335     &                       parcours % gr % Agrif_d(iii))
336            enddo
337C
338            eps = min(g_eps,oldgrid_eps)/100.                 
339C
340            do iii = 1 , Agrif_Probdim
341
342               if (g % Agrif_d(iii) .LT. 
343     &             (parcours % gr % Agrif_d(iii) - eps)) then
344C           
345                   parcours => parcours % next
346C           
347                   out = 1
348C   
349                   Exit
350C             
351               endif
352C     
353            enddo
354        if ( out .EQ. 1 ) Cycle
355C
356            Call AGRIF_CopyFromOld(g,parcours%gr)
357C
358        endif
359C                     
360        Call Agrif_CopyFromold_All
361     &             (g, parcours % gr % child_grids)
362C       
363        parcours => parcours % next
364C       
365      enddo
366C
367C     
368      Return     
369C
370C
371      End Subroutine AGRIF_CopyFromold_All
372C
373C
374C     
375C     **************************************************************************
376CCC   Subroutine Agrif_CopyFromold
377C     **************************************************************************
378C     
379      Subroutine Agrif_CopyFromold(Agrif_New_Gr,Agrif_Old_Gr)
380C
381CCC   Description:
382CCC   Call to the Agrif_Copy procedure.
383C
384CC    Method:
385CC           
386C
387C     Declarations:
388C
389     
390C     
391C     Pointer argument   
392      TYPE(Agrif_Grid),Pointer   :: Agrif_New_Gr  ! Pointer on the current grid
393      TYPE(Agrif_Grid),Pointer :: Agrif_Old_Gr    ! Pointer on an old grid
394      INTEGER :: i
395C
396C     
397      do i = 1 , AGRIF_NbVariables
398         if ( Agrif_Mygrid % tabvars(i) % var % restaure ) then
399C
400            Call Agrif_Copy(Agrif_New_Gr,Agrif_Old_Gr,
401     &           Agrif_New_Gr % tabvars(i),
402     &           Agrif_Old_Gr % tabvars(i))
403C
404        endif
405      enddo
406
407C     
408C           
409      Return
410C
411C     
412      End Subroutine Agrif_CopyFromold
413     
414CC    **************************************************************************
415CCC   Subroutine AGRIF_CopyFromold_AllOneVar
416C     **************************************************************************
417C
418      Recursive Subroutine AGRIF_CopyFromold_AllOneVar(g,oldchildgrids,
419     &   indic)
420C
421CCC   Description:
422CCC   Routine called in the Agrif_Init_Hierarchy procedure 
423C       (Agrif_Clustering module). 
424C
425CC    Method:       
426C
427C     Declarations:
428C
429     
430C     
431C     Pointer argument   
432      TYPE(AGRIF_grid),pointer   :: g ! Pointer on the current grid
433      TYPE(AGRIF_pgrid),pointer   :: oldchildgrids
434      integer :: indic
435C
436C     Local pointer
437      TYPE(AGRIF_pgrid),pointer  :: parcours ! Pointer for the recursive
438                                             ! procedure
439      REAL g_eps,eps,oldgrid_eps
440      INTEGER :: out
441      INTEGER :: iii
442C
443      out = 0
444C                                                                           
445      parcours => oldchildgrids 
446C
447      do while (associated(parcours))
448C 
449        if ((.NOT. g % fixed) .AND. (parcours % gr %oldgrid)) then
450C       
451            g_eps = huge(1.)
452            oldgrid_eps = huge(1.)
453            do iii = 1 , Agrif_Probdim
454               g_eps = min(g_eps,g % Agrif_d(iii))
455               oldgrid_eps = min(oldgrid_eps,
456     &                       parcours % gr % Agrif_d(iii))
457            enddo
458C
459            eps = min(g_eps,oldgrid_eps)/100.                 
460C
461            do iii = 1 , Agrif_Probdim
462
463               if (g % Agrif_d(iii) .LT. 
464     &             (parcours % gr % Agrif_d(iii) - eps)) then
465C           
466                   parcours => parcours % next
467C           
468                   out = 1
469C   
470                   Exit
471C             
472               endif
473C     
474            enddo
475        if ( out .EQ. 1 ) Cycle
476C
477            Call AGRIF_CopyFromOldOneVar(g,parcours%gr,indic)
478C
479        endif
480C                     
481        Call Agrif_CopyFromold_AllOneVar
482     &             (g, parcours % gr % child_grids,indic)
483C       
484        parcours => parcours % next
485C       
486      enddo
487C
488C     
489      Return     
490C
491C
492      End Subroutine AGRIF_CopyFromold_AllOneVar
493C
494C
495C     
496C     **************************************************************************
497CCC   Subroutine Agrif_CopyFromoldOneVar
498C     **************************************************************************
499C     
500      Subroutine Agrif_CopyFromoldOneVar(Agrif_New_Gr,Agrif_Old_Gr,
501     &    indic)
502C
503CCC   Description:
504CCC   Call to the Agrif_Copy procedure.
505C
506CC    Method:
507CC           
508C
509C     Declarations:
510C
511     
512C     
513C     Pointer argument   
514      TYPE(Agrif_Grid),Pointer   :: Agrif_New_Gr  ! Pointer on the current grid
515      TYPE(Agrif_Grid),Pointer :: Agrif_Old_Gr    ! Pointer on an old grid
516      INTEGER :: indic
517      INTEGER :: i
518      TYPE(Agrif_PVariable),Pointer ::tabvars,oldtabvars     
519C
520C     
521      tabvars => Agrif_Search_Variable(Agrif_New_Gr,-indic)
522      oldtabvars => Agrif_Search_Variable(Agrif_Old_Gr,-indic)
523
524      Call Agrif_Nbdim_Allocation(tabvars%var,
525     &  tabvars%var%lb,tabvars%var%ub,
526     &  tabvars%var%nbdim)
527     
528      Call Agrif_Copy(Agrif_New_Gr,Agrif_Old_Gr,
529     &           tabvars,oldtabvars)
530
531
532C     
533C           
534      Return
535C
536C     
537      End Subroutine Agrif_CopyFromoldOneVar     
538     
539C
540C
541CC
542C
543C
544C     **************************************************************************
545CCC   Subroutine Agrif_Copy
546C     **************************************************************************
547C     
548      Subroutine Agrif_Copy(Agrif_New_Gr,Agrif_Old_Gr,New_Var,Old_Var)
549C     
550CCC   Description:
551CCC   Sets arguments of the Agrif_UpdatenD procedures, n being the number of 
552CCC   DIMENSION of the grid variable.
553C
554CC    Method:
555CC           
556C
557C     Declarations:
558C
559     
560C
561C     Pointer argument   
562      TYPE(Agrif_Grid),Pointer   :: Agrif_New_Gr   ! Pointer on the current grid
563      TYPE(Agrif_Grid), Pointer :: Agrif_Old_Gr    ! Pointer on an old grid
564      TYPE(Agrif_Pvariable) :: New_Var, Old_Var
565      INTEGER :: nbdim                  ! Number of dimensions of the current
566                                        ! grid     
567      INTEGER,DIMENSION(6) :: pttabnew  ! Indexes of the first point in the
568                                        ! domain
569      INTEGER,DIMENSION(6) :: petabnew  ! Indexes of the first point in the
570                                        ! domain
571      INTEGER,DIMENSION(6) :: pttabold  ! Indexes of the first point in the
572                                        ! domain
573      INTEGER,DIMENSION(6) :: petabold  ! Indexes of the first point in the
574                                        ! domain
575      INTEGER,DIMENSION(6) :: nbtabold  ! Number of cells in each direction
576     
577      INTEGER,DIMENSION(6) :: nbtabnew  ! Number of cells in each direction   
578      TYPE(AGRIF_Variable), Pointer :: root ! Pointer on the variable of the
579                                            ! root grid
580      REAL, DIMENSION(6) :: snew,sold
581      REAL, DIMENSION(6) :: dsnew,dsold                             
582      REAL :: eps
583      INTEGER :: n
584C
585C 
586      root => New_Var % var % root_var 
587C     
588      nbdim = root % nbdim
589C         
590      do n=1,nbdim
591C 
592        select case(root % interptab(n))
593C       
594        case('x')               
595C
596          pttabnew(n) = root % point(1)
597C       
598          pttabold(n) = root % point(1)
599C       
600          snew(n) = Agrif_New_Gr % Agrif_x(1)
601C
602          sold(n) = Agrif_Old_Gr % Agrif_x(1) 
603C       
604          dsnew(n) = Agrif_New_Gr % Agrif_d(1)
605C
606          dsold(n) = Agrif_Old_Gr % Agrif_d(1)
607C                     
608          if (root % posvar(n) .EQ. 1) then
609C         
610              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(1)
611C         
612              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(1)         
613C       
614            else
615C         
616              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(1) - 1
617C         
618              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(1) - 1
619C
620              snew(n) = snew(n) + dsnew(n)/2.
621C
622              sold(n) = sold(n) + dsold(n)/2.
623C       
624          endif           
625C       
626        case('y')
627C       
628          pttabnew(n) = root % point(2)
629C       
630          pttabold(n) = root % point(2)
631C       
632          snew(n) = Agrif_New_Gr % Agrif_x(2)
633C       
634          sold(n) = Agrif_Old_Gr % Agrif_x(2)
635C       
636          dsnew(n) = Agrif_New_Gr % Agrif_d(2)
637C       
638          dsold(n) = Agrif_Old_Gr % Agrif_d(2)
639C               
640          if (root % posvar(n) .EQ. 1) then
641C       
642              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(2)
643C         
644              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(2)         
645C       
646            else
647C         
648              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(2) - 1
649C         
650              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(2) - 1
651C
652              snew(n) = snew(n) + dsnew(n)/2.
653C
654              sold(n) = sold(n) + dsold(n)/2.                   
655C       
656          endif
657C       
658        case('z')
659C       
660          pttabnew(n) = root % point(3)
661C       
662          pttabold(n) = root % point(3)
663C       
664          snew(n) = Agrif_New_Gr % Agrif_x(3)
665C       
666          sold(n) = Agrif_Old_Gr % Agrif_x(3)
667C       
668          dsnew(n) = Agrif_New_Gr % Agrif_d(3)
669C       
670          dsold(n) = Agrif_Old_Gr % Agrif_d(3)
671C               
672          if (root % posvar(n) .EQ. 1) then
673C       
674              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(3)
675C         
676              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(3)         
677C       
678            else
679C         
680             petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(3) - 1
681C         
682             petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(3) - 1
683C
684             snew(n) = snew(n) + dsnew(n)/2.
685C
686             sold(n) = sold(n) + dsold(n)/2.                   
687C       
688          endif
689C         
690        case('N')
691C       
692          Call Agrif_nbdim_Get_bound(New_Var%var,
693     &                           pttabnew(n),petabnew(n),
694     &                           n,nbdim)
695C       
696          pttabold(n) = pttabnew(n)
697C     
698          petabold(n) = petabnew(n)
699C             
700          snew(n) = 0.
701C     
702          sold(n) = 0. 
703C     
704          dsnew(n) = 1.
705C     
706          dsold(n) = 1.
707C       
708        end select
709C     
710      enddo
711C     
712      do n = 1,nbdim
713C     
714        nbtabnew(n) = petabnew(n) - pttabnew(n)
715C     
716        nbtabold(n) = petabold(n) - pttabold(n)     
717C 
718      enddo     
719C
720      eps = min(minval(dsnew(1:nbdim)),minval(dsold(1:nbdim)))
721C     
722      eps = eps/100.     
723C     
724      do n = 1,nbdim
725C     
726        if (snew(n) .GT. (sold(n)+dsold(n)*nbtabold(n)+eps)) Return
727C     
728        if ((snew(n)+dsnew(n)*nbtabnew(n)-eps) .LT. sold(n)) Return
729C     
730      enddo
731C           
732C
733      Call AGRIF_CopynD
734     &        (New_Var,Old_Var,pttabold,petabold,pttabnew,petabnew,
735     &         sold,snew,dsold,dsnew,nbdim) 
736C     
737C 
738      Return
739C     
740C
741      End Subroutine Agrif_Copy
742C
743C
744C
745C     **************************************************************************
746CCC   Subroutine Agrif_CopynD
747C     **************************************************************************
748C     
749      Subroutine Agrif_CopynD(New_Var,Old_Var,pttabold,petabold,
750     &                        pttabnew,petabnew,sold,snew,dsold,
751     &                        dsnew,nbdim)
752C     
753CCC   Description:
754CCC   Copy of the nD New_Var variable from the nD Old_Var variable.
755C
756CC    Method:
757CC           
758C
759C     Declarations:
760C
761     
762C
763      INTEGER :: nbdim
764      TYPE(Agrif_Pvariable) :: New_Var, Old_Var     
765      INTEGER,DIMENSION(nbdim) :: pttabnew 
766      INTEGER,DIMENSION(nbdim) :: petabnew
767      INTEGER,DIMENSION(nbdim) :: pttabold
768      INTEGER,DIMENSION(nbdim) :: petabold     
769      REAL, DIMENSION(nbdim) :: snew,sold
770      REAL, DIMENSION(nbdim) :: dsnew,dsold 
771C
772      INTEGER :: i,j,k,l,m,n,i0,j0,k0,l0,m0,n0
773C
774      REAL, DIMENSION(nbdim) :: dim_gmin,dim_gmax
775      REAL, DIMENSION(nbdim) :: dim_newmin,dim_newmax
776      REAL, DIMENSION(nbdim) :: dim_min
777      INTEGER, DIMENSION(nbdim) :: ind_gmin,ind_newmin
778      INTEGER, DIMENSION(nbdim) ::  ind_newmax
779C               
780C
781      do i = 1,nbdim
782C     
783        dim_gmin(i) = sold(i)
784        dim_gmax(i) = sold(i) + (petabold(i)-pttabold(i)) * dsold(i)
785C 
786        dim_newmin(i) = snew(i)
787        dim_newmax(i) = snew(i) + (petabnew(i)-pttabnew(i)) * dsnew(i)
788C
789      enddo         
790C
791      do i = 1,nbdim
792C 
793        if (dim_gmax(i) .LT. dim_newmin(i)) Return
794C     
795        if (dim_gmin(i) .GT. dim_newmax(i)) Return
796C
797      enddo
798C
799C
800      do i = 1,nbdim
801C 
802        ind_newmin(i) = pttabnew(i) - floor(-(max(dim_gmin(i),
803     &                      dim_newmin(i))-dim_newmin(i))/dsnew(i)) 
804C     
805        dim_min(i) = snew(i) + (ind_newmin(i)-pttabnew(i))*dsnew(i)
806C     
807        ind_gmin(i) = pttabold(i) + nint((dim_min(i)-
808     &                dim_gmin(i))/dsold(i))     
809C
810        ind_newmax(i) = pttabnew(i) 
811     &                  + int((min(dim_gmax(i),dim_newmax(i))
812     &                         -dim_newmin(i))/dsnew(i))
813C
814      enddo 
815C
816C
817C
818      SELECT CASE (nbdim)
819      CASE (1)
820         i0 = ind_gmin(1)
821         do i = ind_newmin(1),ind_newmax(1)
822            New_Var % var % array1(i) =
823     &      Old_Var % var % array1(i0)
824            New_Var % var % restore1D(i) = 1   
825         i0 = i0 + int(dsnew(1)/dsold(1))
826         enddo
827      CASE (2)
828
829        i0 = ind_gmin(1)
830        do i = ind_newmin(1),ind_newmax(1)
831        j0 = ind_gmin(2)
832        do j = ind_newmin(2),ind_newmax(2)
833           New_Var % var % array2(i,j) =
834     &            Old_Var % var % array2(i0,j0)
835           New_Var % var % restore2D(i,j) = 1
836               j0 = j0 + int(dsnew(2)/dsold(2))
837        enddo
838        i0 = i0 + int(dsnew(1)/dsold(1))
839        enddo
840      CASE (3)
841        i0 = ind_gmin(1)
842        do i = ind_newmin(1),ind_newmax(1)
843        j0 = ind_gmin(2)
844        do j = ind_newmin(2),ind_newmax(2)
845        k0 = ind_gmin(3)
846        do k = ind_newmin(3),ind_newmax(3)
847           New_Var % var % array3(i,j,k) =
848     &                  Old_Var % var % array3(i0,j0,k0)
849           New_Var % var % restore3D(i,j,k) = 1   
850                     k0 = k0 + int(dsnew(3)/dsold(3))
851        enddo
852        j0 = j0 + int(dsnew(2)/dsold(2))
853        enddo
854        i0 = i0 + int(dsnew(1)/dsold(1))
855        enddo
856      CASE (4)
857        i0 = ind_gmin(1)
858        do i = ind_newmin(1),ind_newmax(1)
859        j0 = ind_gmin(2)
860        do j = ind_newmin(2),ind_newmax(2)
861        k0 = ind_gmin(3)
862        do k = ind_newmin(3),ind_newmax(3)
863        l0 = ind_gmin(4)
864        do l = ind_newmin(4),ind_newmax(4) 
865           New_Var % var % array4(i,j,k,l) =
866     &                        Old_Var % var % array4(i0,j0,k0,l0)
867           New_Var % var % restore4D(i,j,k,l) = 1   
868        l0 = l0 + int(dsnew(4)/dsold(4))
869        enddo
870        k0 = k0 + int(dsnew(3)/dsold(3))
871        enddo
872        j0 = j0 + int(dsnew(2)/dsold(2))
873        enddo
874        i0 = i0 + int(dsnew(1)/dsold(1))
875        enddo
876      CASE (5)
877        i0 = ind_gmin(1)
878        do i = ind_newmin(1),ind_newmax(1)
879        j0 = ind_gmin(2)
880        do j = ind_newmin(2),ind_newmax(2)
881        k0 = ind_gmin(3)
882        do k = ind_newmin(3),ind_newmax(3)
883        l0 = ind_gmin(4)
884        do l = ind_newmin(4),ind_newmax(4) 
885        m0 = ind_gmin(5)
886        do m = ind_newmin(5),ind_newmax(5)
887           New_Var % var % array5(i,j,k,l,m) =
888     &                        Old_Var % var % array5(i0,j0,k0,l0,m0)
889           New_Var % var % restore5D(i,j,k,l,m) = 1
890        m0 = m0 + int(dsnew(5)/dsold(5)) 
891        enddo
892        l0 = l0 + int(dsnew(4)/dsold(4))
893        enddo
894        k0 = k0 + int(dsnew(3)/dsold(3))
895        enddo
896        j0 = j0 + int(dsnew(2)/dsold(2))
897        enddo
898        i0 = i0 + int(dsnew(1)/dsold(1))
899        enddo
900      CASE (6)
901        i0 = ind_gmin(1)
902        do i = ind_newmin(1),ind_newmax(1)
903        j0 = ind_gmin(2)
904        do j = ind_newmin(2),ind_newmax(2)
905        k0 = ind_gmin(3)
906        do k = ind_newmin(3),ind_newmax(3)
907        l0 = ind_gmin(4)
908        do l = ind_newmin(4),ind_newmax(4) 
909        m0 = ind_gmin(5)
910        do m = ind_newmin(5),ind_newmax(5)
911        n0 = ind_gmin(6)
912        do n = ind_newmin(6),ind_newmax(6)
913           New_Var % var % array6(i,j,k,l,m,n) =
914     &                        Old_Var % var % array6(i0,j0,k0,l0,m0,n0)
915           New_Var % var % restore6D(i,j,k,l,m,n) = 1
916        n0 = n0 + int(dsnew(6)/dsold(6)) 
917        enddo
918        m0 = m0 + int(dsnew(5)/dsold(5)) 
919        enddo
920        l0 = l0 + int(dsnew(4)/dsold(4))
921        enddo
922        k0 = k0 + int(dsnew(3)/dsold(3))
923        enddo
924        j0 = j0 + int(dsnew(2)/dsold(2))
925        enddo
926        i0 = i0 + int(dsnew(1)/dsold(1))
927        enddo
928      END SELECT
929C           
930      Return
931C
932C
933      End Subroutine Agrif_CopynD                   
934C
935C
936C
937      End module Agrif_Save
Note: See TracBrowser for help on using the repository browser.