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

source: trunk/AGRIF/AGRIF_FILES/modsauv.F @ 907

Last change on this file since 907 was 662, checked in by opalod, 17 years ago

RB: update Agrif internal routines with a new update scheme and performance improvment

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 26.0 KB
Line 
1!
2! $Id$
3!
4C     AGRIF (Adaptive Grid Refinement In Fortran)
5C
6C     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7C                        Christophe Vouland (Christophe.Vouland@imag.fr)   
8C
9C     This program is free software; you can redistribute it and/or modify
10C     it under the terms of the GNU General Public License as published by
11C     the Free Software Foundation; either version 2 of the License, or
12C     (at your option) any later version.
13C
14C     This program is distributed in the hope that it will be useful,
15C     but WITHOUT ANY WARRANTY; without even the implied warranty of
16C     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17C     GNU General Public License for more details.
18C
19C     You should have received a copy of the GNU General Public License
20C     along with this program; if not, write to the Free Software
21C     Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
22C
23C
24C
25CCC   Module Agrif_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
37C
38      IMPLICIT NONE
39C
40      Contains
41C     Definition of procedures contained in this module.
42C
43C
44C
45C     **************************************************************************
46CCC   Subroutine Agrif_Free_data_before
47C     **************************************************************************
48C
49      Subroutine Agrif_Free_data_before(Agrif_Gr)
50C
51CCC   Description:
52CCC   
53C
54CC    Method:
55CC           
56C
57C     Declarations:
58C
59     
60C     
61C     Pointer argument   
62      TYPE(Agrif_Grid),pointer   :: Agrif_Gr ! Pointer on the current grid
63      INTEGER i
64C
65C   
66      do i = 1 , AGRIF_NbVariables
67         if ( .NOT. Agrif_Mygrid % tabvars(i) % var % restaure) then
68C 
69            if (associated(Agrif_Gr%tabvars(i)%var%array1)) then
70               Deallocate(Agrif_Gr%tabvars(i)%var%array1)
71            endif
72            if (associated(Agrif_Gr%tabvars(i)%var%array2)) then
73               Deallocate(Agrif_Gr%tabvars(i)%var%array2)
74            endif
75            if (associated(Agrif_Gr%tabvars(i)%var%array3)) then
76               Deallocate(Agrif_Gr%tabvars(i)%var%array3)
77            endif
78            if (associated(Agrif_Gr%tabvars(i)%var%array4)) then
79               Deallocate(Agrif_Gr%tabvars(i)%var%array4)
80            endif
81            if (associated(Agrif_Gr%tabvars(i)%var%array5)) then
82               Deallocate(Agrif_Gr%tabvars(i)%var%array5)
83            endif
84            if (associated(Agrif_Gr%tabvars(i)%var%array6)) then
85               Deallocate(Agrif_Gr%tabvars(i)%var%array6)
86            endif
87C
88            if (associated(Agrif_Gr%tabvars(i)%var%darray1)) then
89               Deallocate(Agrif_Gr%tabvars(i)%var%darray1)
90            endif
91            if (associated(Agrif_Gr%tabvars(i)%var%darray2)) then
92               Deallocate(Agrif_Gr%tabvars(i)%var%darray2)
93            endif
94            if (associated(Agrif_Gr%tabvars(i)%var%darray3)) then
95               Deallocate(Agrif_Gr%tabvars(i)%var%darray3)
96            endif
97            if (associated(Agrif_Gr%tabvars(i)%var%darray4)) then
98               Deallocate(Agrif_Gr%tabvars(i)%var%darray4)
99            endif
100            if (associated(Agrif_Gr%tabvars(i)%var%darray5)) then
101               Deallocate(Agrif_Gr%tabvars(i)%var%darray5)
102            endif
103            if (associated(Agrif_Gr%tabvars(i)%var%darray6)) then
104               Deallocate(Agrif_Gr%tabvars(i)%var%darray6)
105            endif
106C
107            if (associated(Agrif_Gr%tabvars(i)%var%larray1)) then
108               Deallocate(Agrif_Gr%tabvars(i)%var%larray1)
109            endif
110            if (associated(Agrif_Gr%tabvars(i)%var%larray2)) then
111               Deallocate(Agrif_Gr%tabvars(i)%var%larray2)
112            endif
113            if (associated(Agrif_Gr%tabvars(i)%var%larray3)) then
114               Deallocate(Agrif_Gr%tabvars(i)%var%larray3)
115            endif
116            if (associated(Agrif_Gr%tabvars(i)%var%larray4)) then
117               Deallocate(Agrif_Gr%tabvars(i)%var%larray4)
118            endif
119            if (associated(Agrif_Gr%tabvars(i)%var%larray5)) then
120               Deallocate(Agrif_Gr%tabvars(i)%var%larray5)
121            endif
122            if (associated(Agrif_Gr%tabvars(i)%var%larray6)) then
123               Deallocate(Agrif_Gr%tabvars(i)%var%larray6)
124            endif
125C
126            if (associated(Agrif_Gr%tabvars(i)%var%iarray1)) then
127               Deallocate(Agrif_Gr%tabvars(i)%var%iarray1)
128            endif
129            if (associated(Agrif_Gr%tabvars(i)%var%iarray2)) then
130               Deallocate(Agrif_Gr%tabvars(i)%var%iarray2)
131            endif
132            if (associated(Agrif_Gr%tabvars(i)%var%iarray3)) then
133               Deallocate(Agrif_Gr%tabvars(i)%var%iarray3)
134            endif
135            if (associated(Agrif_Gr%tabvars(i)%var%iarray4)) then
136               Deallocate(Agrif_Gr%tabvars(i)%var%iarray4)
137            endif
138            if (associated(Agrif_Gr%tabvars(i)%var%iarray5)) then
139               Deallocate(Agrif_Gr%tabvars(i)%var%iarray5)
140            endif
141            if (associated(Agrif_Gr%tabvars(i)%var%iarray6)) then
142               Deallocate(Agrif_Gr%tabvars(i)%var%iarray6)
143            endif
144C
145            if (associated(Agrif_Gr%tabvars(i)%var%carray1)) then
146               Deallocate(Agrif_Gr%tabvars(i)%var%carray1)
147            endif
148            if (associated(Agrif_Gr%tabvars(i)%var%carray2)) then
149               Deallocate(Agrif_Gr%tabvars(i)%var%carray2)
150            endif
151C
152            if (associated(Agrif_Gr%tabvars(i)%var%oldvalues2D)) then
153               Deallocate(Agrif_Gr%tabvars(i)%var%oldvalues2D)
154            endif
155            if (associated(Agrif_Gr%tabvars(i)%var%interpIndex)) then
156               Deallocate(Agrif_Gr%tabvars(i)%var%interpIndex)
157            endif
158           
159            if (associated(Agrif_Gr%tabvars(i)%var%posvar)) then
160               Deallocate(Agrif_Gr%tabvars(i)%var%posvar)
161            endif     
162           
163            if (associated(Agrif_Gr%tabvars(i)%var%interptab)) then
164               Deallocate(Agrif_Gr%tabvars(i)%var%interptab)
165            endif 
166                             
167C
168            Deallocate(Agrif_Gr%tabvars(i)%var)
169C 
170        endif
171      enddo
172C
173C 
174C
175      if ( Agrif_USE_ONLY_FIXED_GRIDS .EQ. 0 ) then
176         if ( Agrif_Probdim .EQ. 1 ) Deallocate(Agrif_Gr%tabpoint1D)
177         if ( Agrif_Probdim .EQ. 2 ) Deallocate(Agrif_Gr%tabpoint2D)
178         if ( Agrif_Probdim .EQ. 3 ) Deallocate(Agrif_Gr%tabpoint3D)
179      endif
180C         
181      Return
182C
183C
184      End Subroutine Agrif_Free_data_before
185C
186C
187C
188C     **************************************************************************
189CCC   Subroutine Agrif_Free_data_after
190C     **************************************************************************
191C     
192      Subroutine Agrif_Free_data_after(Agrif_Gr)
193C
194CCC   Description:
195CCC   
196C
197CC    Method:
198CC           
199C
200C     Declarations:
201C
202C
203     
204C
205C     Pointer argument   
206      TYPE(Agrif_Grid),pointer   :: Agrif_Gr  ! Pointer on the current grid
207      INTEGER i
208C
209C     
210      do i = 1 , AGRIF_NbVariables
211         if ( Agrif_Mygrid % tabvars(i) % var % restaure) then
212C 
213            if (associated(Agrif_Gr%tabvars(i)%var%array1)) then
214               Deallocate(Agrif_Gr%tabvars(i)%var%array1)
215            endif
216            if (associated(Agrif_Gr%tabvars(i)%var%array2)) then
217               Deallocate(Agrif_Gr%tabvars(i)%var%array2)
218            endif
219            if (associated(Agrif_Gr%tabvars(i)%var%array3)) then
220               Deallocate(Agrif_Gr%tabvars(i)%var%array3)
221            endif
222            if (associated(Agrif_Gr%tabvars(i)%var%array4)) then
223               Deallocate(Agrif_Gr%tabvars(i)%var%array4)
224            endif
225            if (associated(Agrif_Gr%tabvars(i)%var%array5)) then
226               Deallocate(Agrif_Gr%tabvars(i)%var%array5)
227            endif
228            if (associated(Agrif_Gr%tabvars(i)%var%array6)) then
229               Deallocate(Agrif_Gr%tabvars(i)%var%array6)
230            endif
231!
232            if (associated(Agrif_Gr%tabvars(i)%var%darray1)) then
233               Deallocate(Agrif_Gr%tabvars(i)%var%darray1)
234            endif
235            if (associated(Agrif_Gr%tabvars(i)%var%darray2)) then
236               Deallocate(Agrif_Gr%tabvars(i)%var%darray2)
237            endif
238            if (associated(Agrif_Gr%tabvars(i)%var%darray3)) then
239               Deallocate(Agrif_Gr%tabvars(i)%var%darray3)
240            endif
241            if (associated(Agrif_Gr%tabvars(i)%var%darray4)) then
242               Deallocate(Agrif_Gr%tabvars(i)%var%darray4)
243            endif
244            if (associated(Agrif_Gr%tabvars(i)%var%darray5)) then
245               Deallocate(Agrif_Gr%tabvars(i)%var%darray5)
246            endif
247            if (associated(Agrif_Gr%tabvars(i)%var%darray6)) then
248               Deallocate(Agrif_Gr%tabvars(i)%var%darray6)
249            endif
250!
251            if (associated(Agrif_Gr%tabvars(i)%var%larray1)) then
252               Deallocate(Agrif_Gr%tabvars(i)%var%larray1)
253            endif
254            if (associated(Agrif_Gr%tabvars(i)%var%larray2)) then
255               Deallocate(Agrif_Gr%tabvars(i)%var%larray2)
256            endif
257            if (associated(Agrif_Gr%tabvars(i)%var%larray3)) then
258               Deallocate(Agrif_Gr%tabvars(i)%var%larray3)
259            endif
260            if (associated(Agrif_Gr%tabvars(i)%var%larray4)) then
261               Deallocate(Agrif_Gr%tabvars(i)%var%larray4)
262            endif
263            if (associated(Agrif_Gr%tabvars(i)%var%larray5)) then
264               Deallocate(Agrif_Gr%tabvars(i)%var%larray5)
265            endif
266            if (associated(Agrif_Gr%tabvars(i)%var%larray6)) then
267               Deallocate(Agrif_Gr%tabvars(i)%var%larray6)
268            endif
269!
270            if (associated(Agrif_Gr%tabvars(i)%var%iarray1)) then
271               Deallocate(Agrif_Gr%tabvars(i)%var%iarray1)
272            endif
273            if (associated(Agrif_Gr%tabvars(i)%var%iarray2)) then
274               Deallocate(Agrif_Gr%tabvars(i)%var%iarray2)
275            endif
276            if (associated(Agrif_Gr%tabvars(i)%var%iarray3)) then
277               Deallocate(Agrif_Gr%tabvars(i)%var%iarray3)
278            endif
279            if (associated(Agrif_Gr%tabvars(i)%var%iarray4)) then
280               Deallocate(Agrif_Gr%tabvars(i)%var%iarray4)
281            endif
282            if (associated(Agrif_Gr%tabvars(i)%var%iarray5)) then
283               Deallocate(Agrif_Gr%tabvars(i)%var%iarray5)
284            endif
285            if (associated(Agrif_Gr%tabvars(i)%var%iarray6)) then
286               Deallocate(Agrif_Gr%tabvars(i)%var%iarray6)
287            endif
288!
289            if (associated(Agrif_Gr%tabvars(i)%var%carray1)) then
290               Deallocate(Agrif_Gr%tabvars(i)%var%carray1)
291            endif
292            if (associated(Agrif_Gr%tabvars(i)%var%carray2)) then
293               Deallocate(Agrif_Gr%tabvars(i)%var%carray2)
294            endif
295!
296            if (associated(Agrif_Gr%tabvars(i)%var%oldvalues2D)) then
297               Deallocate(Agrif_Gr%tabvars(i)%var%oldvalues2D)
298            endif
299            if (associated(Agrif_Gr%tabvars(i)%var%interpIndex)) then
300               Deallocate(Agrif_Gr%tabvars(i)%var%interpIndex)
301            endif
302           
303            if (associated(Agrif_Gr%tabvars(i)%var%posvar)) then
304               Deallocate(Agrif_Gr%tabvars(i)%var%posvar)
305            endif     
306           
307            if (associated(Agrif_Gr%tabvars(i)%var%interptab)) then
308               Deallocate(Agrif_Gr%tabvars(i)%var%interptab)
309            endif             
310!
311            Deallocate(Agrif_Gr%tabvars(i)%var)
312!
313         endif
314      enddo
315C
316C     
317      Return
318C
319C
320      End Subroutine Agrif_Free_data_after     
321C
322C
323CC    **************************************************************************
324CCC   Subroutine AGRIF_CopyFromold_All
325C     **************************************************************************
326C
327      Recursive Subroutine AGRIF_CopyFromold_All(g,oldchildgrids)
328C
329CCC   Description:
330CCC   Routine called in the Agrif_Init_Hierarchy procedure 
331C       (Agrif_Clustering module). 
332C
333CC    Method:       
334C
335C     Declarations:
336C
337     
338C     
339C     Pointer argument   
340      TYPE(AGRIF_grid),pointer   :: g ! Pointer on the current grid
341      TYPE(AGRIF_pgrid),pointer   :: oldchildgrids
342C
343C     Local pointer
344      TYPE(AGRIF_pgrid),pointer  :: parcours ! Pointer for the recursive
345                                             ! procedure
346      REAL g_eps,eps,oldgrid_eps
347      INTEGER :: out
348      INTEGER :: iii
349C
350      out = 0
351C                                                                           
352      parcours => oldchildgrids 
353C
354      do while (associated(parcours))
355C 
356        if ((.NOT. g % fixed) .AND. (parcours % gr %oldgrid)) then
357C       
358            g_eps = huge(1.)
359            oldgrid_eps = huge(1.)
360            do iii = 1 , Agrif_Probdim
361               g_eps = min(g_eps,g % Agrif_d(iii))
362               oldgrid_eps = min(oldgrid_eps,
363     &                       parcours % gr % Agrif_d(iii))
364            enddo
365C
366            eps = min(g_eps,oldgrid_eps)/100.                 
367C
368            do iii = 1 , Agrif_Probdim
369
370               if (g % Agrif_d(iii) .LT. 
371     &             (parcours % gr % Agrif_d(iii) - eps)) then
372C           
373                   parcours => parcours % next
374C           
375                   out = 1
376C   
377                   Exit
378C             
379               endif
380C     
381            enddo
382        if ( out .EQ. 1 ) Cycle
383C
384            Call AGRIF_CopyFromOld(g,parcours%gr)
385C
386        endif
387C                     
388        Call Agrif_CopyFromold_All
389     &             (g, parcours % gr % child_grids)
390C       
391        parcours => parcours % next
392C       
393      enddo
394C
395C     
396      Return     
397C
398C
399      End Subroutine AGRIF_CopyFromold_All
400C
401C
402C     
403C     **************************************************************************
404CCC   Subroutine Agrif_CopyFromold
405C     **************************************************************************
406C     
407      Subroutine Agrif_CopyFromold(Agrif_New_Gr,Agrif_Old_Gr)
408C
409CCC   Description:
410CCC   Call to the Agrif_Copy procedure.
411C
412CC    Method:
413CC           
414C
415C     Declarations:
416C
417     
418C     
419C     Pointer argument   
420      TYPE(Agrif_Grid),Pointer   :: Agrif_New_Gr  ! Pointer on the current grid
421      TYPE(Agrif_Grid),Pointer :: Agrif_Old_Gr    ! Pointer on an old grid
422      INTEGER :: i
423C
424C     
425      do i = 1 , AGRIF_NbVariables
426         if ( Agrif_Mygrid % tabvars(i) % var % restaure ) then
427C
428            Call Agrif_Copy(Agrif_New_Gr,Agrif_Old_Gr,
429     &           Agrif_New_Gr % tabvars(i),
430     &           Agrif_Old_Gr % tabvars(i))
431C
432        endif
433      enddo
434
435C     
436C           
437      Return
438C
439C     
440      End Subroutine Agrif_CopyFromold
441C
442C
443CC
444C
445C
446C     **************************************************************************
447CCC   Subroutine Agrif_Copy
448C     **************************************************************************
449C     
450      Subroutine Agrif_Copy(Agrif_New_Gr,Agrif_Old_Gr,New_Var,Old_Var)
451C     
452CCC   Description:
453CCC   Sets arguments of the Agrif_UpdatenD procedures, n being the number of 
454CCC   DIMENSION of the grid variable.
455C
456CC    Method:
457CC           
458C
459C     Declarations:
460C
461     
462C
463C     Pointer argument   
464      TYPE(Agrif_Grid),Pointer   :: Agrif_New_Gr   ! Pointer on the current grid
465      TYPE(Agrif_Grid), Pointer :: Agrif_Old_Gr    ! Pointer on an old grid
466      TYPE(Agrif_Pvariable) :: New_Var, Old_Var
467      INTEGER :: nbdim                  ! Number of dimensions of the current
468                                        ! grid     
469      INTEGER,DIMENSION(6) :: pttabnew  ! Indexes of the first point in the
470                                        ! domain
471      INTEGER,DIMENSION(6) :: petabnew  ! Indexes of the first point in the
472                                        ! domain
473      INTEGER,DIMENSION(6) :: pttabold  ! Indexes of the first point in the
474                                        ! domain
475      INTEGER,DIMENSION(6) :: petabold  ! Indexes of the first point in the
476                                        ! domain
477      INTEGER,DIMENSION(6) :: nbtabold  ! Number of cells in each direction
478     
479      INTEGER,DIMENSION(6) :: nbtabnew  ! Number of cells in each direction   
480      TYPE(AGRIF_Variable), Pointer :: root ! Pointer on the variable of the
481                                            ! root grid
482      REAL, DIMENSION(6) :: snew,sold
483      REAL, DIMENSION(6) :: dsnew,dsold                             
484      REAL :: eps
485      INTEGER :: n
486C
487C 
488      root => New_Var % var % root_var 
489C     
490      nbdim = root % nbdim
491C         
492      do n=1,nbdim
493C 
494        select case(root % interptab(n))
495C       
496        case('x')               
497C
498          pttabnew(n) = root % point(1)
499C       
500          pttabold(n) = root % point(1)
501C       
502          snew(n) = Agrif_New_Gr % Agrif_x(1)
503C
504          sold(n) = Agrif_Old_Gr % Agrif_x(1) 
505C       
506          dsnew(n) = Agrif_New_Gr % Agrif_d(1)
507C
508          dsold(n) = Agrif_Old_Gr % Agrif_d(1)
509C                     
510          if (root % posvar(n) .EQ. 1) then
511C         
512              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(1)
513C         
514              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(1)         
515C       
516            else
517C         
518              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(1) - 1
519C         
520              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(1) - 1
521C
522              snew(n) = snew(n) + dsnew(n)/2.
523C
524              sold(n) = sold(n) + dsold(n)/2.
525C       
526          endif           
527C       
528        case('y')
529C       
530          pttabnew(n) = root % point(2)
531C       
532          pttabold(n) = root % point(2)
533C       
534          snew(n) = Agrif_New_Gr % Agrif_x(2)
535C       
536          sold(n) = Agrif_Old_Gr % Agrif_x(2)
537C       
538          dsnew(n) = Agrif_New_Gr % Agrif_d(2)
539C       
540          dsold(n) = Agrif_Old_Gr % Agrif_d(2)
541C               
542          if (root % posvar(n) .EQ. 1) then
543C       
544              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(2)
545C         
546              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(2)         
547C       
548            else
549C         
550              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(2) - 1
551C         
552              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(2) - 1
553C
554              snew(n) = snew(n) + dsnew(n)/2.
555C
556              sold(n) = sold(n) + dsold(n)/2.                   
557C       
558          endif
559C       
560        case('z')
561C       
562          pttabnew(n) = root % point(3)
563C       
564          pttabold(n) = root % point(3)
565C       
566          snew(n) = Agrif_New_Gr % Agrif_x(3)
567C       
568          sold(n) = Agrif_Old_Gr % Agrif_x(3)
569C       
570          dsnew(n) = Agrif_New_Gr % Agrif_d(3)
571C       
572          dsold(n) = Agrif_Old_Gr % Agrif_d(3)
573C               
574          if (root % posvar(n) .EQ. 1) then
575C       
576              petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(3)
577C         
578              petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(3)         
579C       
580            else
581C         
582             petabnew(n) = pttabnew(n) + Agrif_New_Gr % nb(3) - 1
583C         
584             petabold(n) = pttabold(n) + Agrif_Old_Gr % nb(3) - 1
585C
586             snew(n) = snew(n) + dsnew(n)/2.
587C
588             sold(n) = sold(n) + dsold(n)/2.                   
589C       
590          endif
591C         
592        case('N')
593C       
594          Call Agrif_nbdim_Get_bound(New_Var%var,
595     &                           pttabnew(n),petabnew(n),
596     &                           n,nbdim)
597C       
598          pttabold(n) = pttabnew(n)
599C     
600          petabold(n) = petabnew(n)
601C             
602          snew(n) = 0.
603C     
604          sold(n) = 0. 
605C     
606          dsnew(n) = 1.
607C     
608          dsold(n) = 1.
609C       
610        end select
611C     
612      enddo
613C     
614      do n = 1,nbdim
615C     
616        nbtabnew(n) = petabnew(n) - pttabnew(n)
617C     
618        nbtabold(n) = petabold(n) - pttabold(n)     
619C 
620      enddo     
621C
622      eps = min(minval(dsnew(1:nbdim)),minval(dsold(1:nbdim)))
623C     
624      eps = eps/100.     
625C     
626      do n = 1,nbdim
627C     
628        if (snew(n) .GT. (sold(n)+dsold(n)*nbtabold(n)+eps)) Return
629C     
630        if ((snew(n)+dsnew(n)*nbtabnew(n)-eps) .LT. sold(n)) Return
631C     
632      enddo
633C           
634C
635      Call AGRIF_CopynD
636     &        (New_Var,Old_Var,pttabold,petabold,pttabnew,petabnew,
637     &         sold,snew,dsold,dsnew,nbdim) 
638C     
639C 
640      Return
641C     
642C
643      End Subroutine Agrif_Copy
644C
645C
646C
647C     **************************************************************************
648CCC   Subroutine Agrif_CopynD
649C     **************************************************************************
650C     
651      Subroutine Agrif_CopynD(New_Var,Old_Var,pttabold,petabold,
652     &                        pttabnew,petabnew,sold,snew,dsold,
653     &                        dsnew,nbdim)
654C     
655CCC   Description:
656CCC   Copy of the nD New_Var variable from the nD Old_Var variable.
657C
658CC    Method:
659CC           
660C
661C     Declarations:
662C
663     
664C
665      INTEGER :: nbdim
666      TYPE(Agrif_Pvariable) :: New_Var, Old_Var     
667      INTEGER,DIMENSION(nbdim) :: pttabnew 
668      INTEGER,DIMENSION(nbdim) :: petabnew
669      INTEGER,DIMENSION(nbdim) :: pttabold
670      INTEGER,DIMENSION(nbdim) :: petabold     
671      REAL, DIMENSION(nbdim) :: snew,sold
672      REAL, DIMENSION(nbdim) :: dsnew,dsold 
673C
674      INTEGER :: i,j,k,l,m,n,i0,j0,k0,l0,m0,n0
675C
676      REAL, DIMENSION(nbdim) :: dim_gmin,dim_gmax
677      REAL, DIMENSION(nbdim) :: dim_newmin,dim_newmax
678      REAL, DIMENSION(nbdim) :: dim_min
679      INTEGER, DIMENSION(nbdim) :: ind_gmin,ind_newmin
680      INTEGER, DIMENSION(nbdim) ::  ind_newmax
681C               
682C
683      do i = 1,nbdim
684C     
685        dim_gmin(i) = sold(i)
686        dim_gmax(i) = sold(i) + (petabold(i)-pttabold(i)) * dsold(i)
687C 
688        dim_newmin(i) = snew(i)
689        dim_newmax(i) = snew(i) + (petabnew(i)-pttabnew(i)) * dsnew(i)
690C
691      enddo         
692C
693      do i = 1,nbdim
694C 
695        if (dim_gmax(i) .LT. dim_newmin(i)) Return
696C     
697        if (dim_gmin(i) .GT. dim_newmax(i)) Return
698C
699      enddo
700C
701C
702      do i = 1,nbdim
703C 
704        ind_newmin(i) = pttabnew(i) - floor(-(max(dim_gmin(i),
705     &                      dim_newmin(i))-dim_newmin(i))/dsnew(i)) 
706C     
707        dim_min(i) = snew(i) + (ind_newmin(i)-pttabnew(i))*dsnew(i)
708C     
709        ind_gmin(i) = pttabold(i) + nint((dim_min(i)-
710     &                dim_gmin(i))/dsold(i))     
711C
712        ind_newmax(i) = pttabnew(i) 
713     &                  + int((min(dim_gmax(i),dim_newmax(i))
714     &                         -dim_newmin(i))/dsnew(i))
715C
716      enddo 
717C
718C
719C
720      SELECT CASE (nbdim)
721      CASE (1)
722         i0 = ind_gmin(1)
723         do i = ind_newmin(1),ind_newmax(1)
724            New_Var % var % array1(i) =
725     &      Old_Var % var % array1(i0)
726            New_Var % var % restore1D(i) = 1   
727         i0 = i0 + int(dsnew(1)/dsold(1))
728         enddo
729      CASE (2)
730        i0 = ind_gmin(1)
731        do i = ind_newmin(1),ind_newmax(1)
732        j0 = ind_gmin(2)
733        do j = ind_newmin(2),ind_newmax(2)
734           New_Var % var % array2(i,j) =
735     &            Old_Var % var % array2(i0,j0)
736           New_Var % var % restore2D(i,j) = 1
737               j0 = j0 + int(dsnew(2)/dsold(2))
738        enddo
739        i0 = i0 + int(dsnew(1)/dsold(1))
740        enddo
741      CASE (3)
742        i0 = ind_gmin(1)
743        do i = ind_newmin(1),ind_newmax(1)
744        j0 = ind_gmin(2)
745        do j = ind_newmin(2),ind_newmax(2)
746        k0 = ind_gmin(3)
747        do k = ind_newmin(3),ind_newmax(3)
748           New_Var % var % array3(i,j,k) =
749     &                  Old_Var % var % array3(i0,j0,k0)
750           New_Var % var % restore3D(i,j,k) = 1   
751                     k0 = k0 + int(dsnew(3)/dsold(3))
752        enddo
753        j0 = j0 + int(dsnew(2)/dsold(2))
754        enddo
755        i0 = i0 + int(dsnew(1)/dsold(1))
756        enddo
757      CASE (4)
758        i0 = ind_gmin(1)
759        do i = ind_newmin(1),ind_newmax(1)
760        j0 = ind_gmin(2)
761        do j = ind_newmin(2),ind_newmax(2)
762        k0 = ind_gmin(3)
763        do k = ind_newmin(3),ind_newmax(3)
764        l0 = ind_gmin(4)
765        do l = ind_newmin(4),ind_newmax(4) 
766           New_Var % var % array4(i,j,k,l) =
767     &                        Old_Var % var % array4(i0,j0,k0,l0)
768           New_Var % var % restore4D(i,j,k,l) = 1   
769        l0 = l0 + int(dsnew(4)/dsold(4))
770        enddo
771        k0 = k0 + int(dsnew(3)/dsold(3))
772        enddo
773        j0 = j0 + int(dsnew(2)/dsold(2))
774        enddo
775        i0 = i0 + int(dsnew(1)/dsold(1))
776        enddo
777      CASE (5)
778        i0 = ind_gmin(1)
779        do i = ind_newmin(1),ind_newmax(1)
780        j0 = ind_gmin(2)
781        do j = ind_newmin(2),ind_newmax(2)
782        k0 = ind_gmin(3)
783        do k = ind_newmin(3),ind_newmax(3)
784        l0 = ind_gmin(4)
785        do l = ind_newmin(4),ind_newmax(4) 
786        m0 = ind_gmin(5)
787        do m = ind_newmin(5),ind_newmax(5)
788           New_Var % var % array5(i,j,k,l,m) =
789     &                        Old_Var % var % array5(i0,j0,k0,l0,m0)
790           New_Var % var % restore5D(i,j,k,l,m) = 1
791        m0 = m0 + int(dsnew(5)/dsold(5)) 
792        enddo
793        l0 = l0 + int(dsnew(4)/dsold(4))
794        enddo
795        k0 = k0 + int(dsnew(3)/dsold(3))
796        enddo
797        j0 = j0 + int(dsnew(2)/dsold(2))
798        enddo
799        i0 = i0 + int(dsnew(1)/dsold(1))
800        enddo
801      CASE (6)
802        i0 = ind_gmin(1)
803        do i = ind_newmin(1),ind_newmax(1)
804        j0 = ind_gmin(2)
805        do j = ind_newmin(2),ind_newmax(2)
806        k0 = ind_gmin(3)
807        do k = ind_newmin(3),ind_newmax(3)
808        l0 = ind_gmin(4)
809        do l = ind_newmin(4),ind_newmax(4) 
810        m0 = ind_gmin(5)
811        do m = ind_newmin(5),ind_newmax(5)
812        n0 = ind_gmin(6)
813        do n = ind_newmin(6),ind_newmax(6)
814           New_Var % var % array6(i,j,k,l,m,n) =
815     &                        Old_Var % var % array6(i0,j0,k0,l0,m0,n0)
816           New_Var % var % restore6D(i,j,k,l,m,n) = 1
817        n0 = n0 + int(dsnew(6)/dsold(6)) 
818        enddo
819        m0 = m0 + int(dsnew(5)/dsold(5)) 
820        enddo
821        l0 = l0 + int(dsnew(4)/dsold(4))
822        enddo
823        k0 = k0 + int(dsnew(3)/dsold(3))
824        enddo
825        j0 = j0 + int(dsnew(2)/dsold(2))
826        enddo
827        i0 = i0 + int(dsnew(1)/dsold(1))
828        enddo
829      END SELECT
830C           
831      Return
832C
833C
834      End Subroutine Agrif_CopynD                   
835C
836C
837C
838      End module Agrif_Save
Note: See TracBrowser for help on using the repository browser.