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

source: trunk/AGRIF/AGRIF_FILES/modbc.F @ 662

Last change on this file since 662 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: 46.7 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_Boundary 
26C 
27      Module Agrif_Boundary 
28C
29CCC   Description:
30CCC   Module to calculate the boundary conditions on the child grids from their
31CCC   parent grids.
32C
33C     Modules used: 
34C   
35      Use Agrif_Interpolation       
36C
37      IMPLICIT NONE
38C       
39      CONTAINS   
40C     Define procedures contained in this module
41C
42C
43C
44C     **************************************************************************
45CCC   Subroutine Agrif_Interp_bc_1d
46C     **************************************************************************
47C 
48      Subroutine Agrif_Interp_bc_1d(TypeInterp,parent,child,tab,deb,fin,
49     &                              weight,pweight)           
50C
51CCC   Description:
52CCC   Subroutine to calculate the boundary conditions on a fine grid for a 1D 
53CCC   grid variable.
54C
55C     Declarations:
56C     
57     
58C
59C     Arguments     
60      INTEGER,DIMENSION(6) :: TypeInterp    ! TYPE of interpolation
61                                            ! (linear,...)
62      TYPE(AGRIF_PVariable) :: parent       ! Variable on the parent grid
63      TYPE(AGRIF_PVariable) :: child        ! Variable on the child grid
64      TYPE(AGRIF_PVariable) :: childtemp    ! Temporary variable on the child
65                                            ! grid
66      INTEGER :: deb,fin                    ! Positions of the  interpolations
67      REAL, DIMENSION(
68     &    lbound(child%var%array1,1):ubound(child%var%array1,1)
69     &    ), Target :: tab ! Values of the grid variable
70      LOGICAL :: pweight                    ! Indicates if weight is used for 
71                                            ! the temporal interpolation
72      REAL :: weight                        ! Coefficient for the time
73                                            ! interpolation
74C
75C
76C     Definition of a temporary AGRIF_PVariable data TYPE representing the grid
77C     variable. 
78C
79      allocate(childtemp % var) 
80C
81      childtemp % var % root_var => child % var % root_var
82C     
83C     Values of the grid variable
84      childtemp % var % array1 => tab 
85C
86C     Temporary results for the time interpolation before and after the space 
87C     interpolation 
88      childtemp % var % oldvalues2D => child % var % oldvalues2D
89C 
90C     Index indicating if a space interpolation is necessary
91      childtemp % var % interpIndex => child % var % interpIndex       
92      childtemp % var % Interpolationshouldbemade = 
93     &                 child % var % Interpolationshouldbemade 
94      childtemp % var % list_interp => child % var% list_interp         
95C
96C     Call to the procedure for the calculations of the boundary conditions
97      Call Agrif_CorrectVariable
98     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight)
99C
100      child % var % oldvalues2D => childtemp % var % oldvalues2D
101      child % var % list_interp => childtemp % var %list_interp     
102C     
103      deallocate(childtemp % var)
104C
105C       
106      End Subroutine Agrif_Interp_bc_1D
107C
108C
109C
110C     **************************************************************************
111CCC   Subroutine Agrif_Interp_bc_2d
112C     **************************************************************************
113C 
114      Subroutine Agrif_Interp_bc_2d(TypeInterp,parent,child,tab,deb,fin,
115     &                              weight,pweight,procname)
116C
117CCC   Description:
118CCC   Subroutine to calculate the boundary conditions on a fine grid for a 2D 
119CCC   grid variable.
120C
121C     Declarations:
122C     
123     
124C
125C     Arguments     
126      External :: procname
127      Optional :: procname
128      INTEGER,DIMENSION(6) :: TypeInterp    ! TYPE of interpolation (linear,
129                                            ! lagrange, spline, ... )
130      TYPE(AGRIF_PVariable) :: parent       ! Variable on the parent grid
131      TYPE(AGRIF_PVariable) :: child        ! Variable on the child grid
132      TYPE(AGRIF_PVariable) :: childtemp    ! Temporary variable on the child
133                                            ! grid
134      INTEGER :: deb,fin                    ! Positions where interpolations are
135                                            ! done on the fine grid
136      REAL, DIMENSION(
137     &         lbound(child%var%array2,1): ubound(child%var%array2,1),
138     &         lbound(child%var%array2,2): ubound(child%var%array2,2)), 
139     &         Target :: tab ! Values of the grid variable
140      LOGICAL :: pweight                    ! Indicates if weight is used for 
141                                            ! the temporal interpolation
142      REAL :: weight                        ! Coefficient for the time
143                                            ! interpolation
144C
145C
146C     Definition of a temporary AGRIF_PVariable data TYPE representing the grid
147C     variable.   
148C
149      allocate(childtemp % var) 
150C
151      childtemp % var % root_var => child % var % root_var
152C     
153C     Values of the grid variable
154      childtemp % var % array2 => tab 
155C
156C     Temporary results for the time interpolation before and after the space 
157C     interpolation 
158      childtemp % var % oldvalues2D => child % var % oldvalues2D
159C 
160C     Index indicating if a space interpolation is necessary
161      childtemp % var % interpIndex => child % var % interpIndex       
162      childtemp % var % Interpolationshouldbemade = 
163     &                 child % var % Interpolationshouldbemade   
164      childtemp % var % list_interp => child % var% list_interp   
165C
166C     Call to the procedure for the calculations of the boundary conditions
167      IF (present(procname)) THEN
168      Call Agrif_CorrectVariable
169     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname)
170      ELSE
171      Call Agrif_CorrectVariable
172     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight)
173      ENDIF
174
175C   
176      child % var % oldvalues2D => childtemp % var % oldvalues2D
177      child % var % list_interp => childtemp % var %list_interp
178C         
179      deallocate(childtemp % var)
180C
181C       
182      End Subroutine Agrif_Interp_bc_2D
183C
184C
185C
186C     **************************************************************************
187CCC   Subroutine Agrif_Interp_bc_3d
188C     **************************************************************************
189C 
190      Subroutine Agrif_Interp_bc_3d(TypeInterp,parent,child,tab,deb,fin,
191     &                              weight,pweight,procname)           
192C
193CCC   Description:
194CCC   Subroutine to calculate the boundary conditions on a fine grid for a 3D 
195CCC   variable.
196C
197C     Declarations:
198C     
199     
200C
201C     Arguments 
202      External :: procname
203      Optional :: procname     
204      INTEGER,DIMENSION(6) :: TypeInterp    ! TYPE of interpolation (linear,
205                                            ! lagrange, spline, ... )
206      TYPE(AGRIF_PVariable) :: parent       ! Variable on the parent grid
207      TYPE(AGRIF_PVariable) :: child        ! Variable on the child grid
208      TYPE(AGRIF_PVariable) :: childtemp    ! Temporary variable on the child
209                                            ! grid
210      INTEGER :: deb,fin                    ! Positions where interpolations
211                                            ! are done on the fine grid
212      REAL, DIMENSION(
213     &         lbound(child%var%array3,1):ubound(child%var%array3,1),
214     &         lbound(child%var%array3,2):ubound(child%var%array3,2),
215     &         lbound(child%var%array3,3):ubound(child%var%array3,3)
216     &         ), Target :: tab ! Values of the grid variable
217      LOGICAL :: pweight                    ! Indicates if weight is used for
218                                            ! the temporal interpolation
219      REAL :: weight                        ! Coefficient for the time
220                                            ! interpolation
221C
222C
223C     Definition of a temporary AGRIF_PVariable data TYPE representing the grid 
224C     variable. 
225C
226      allocate(childtemp % var) 
227C
228      childtemp % var % root_var => child % var % root_var       
229C     
230C     Values of the grid variable
231      childtemp % var % array3 => tab
232C
233C     Temporary results for the time interpolation before and after the space 
234C     interpolation 
235      childtemp % var % oldvalues2D => child % var % oldvalues2D
236C 
237C     Index indicating if a space interpolation is necessary
238      childtemp % var % interpIndex => child % var % interpIndex
239      childtemp % var % Interpolationshouldbemade = 
240     &                 child % var % Interpolationshouldbemade 
241      childtemp % var % list_interp => child % var% list_interp           
242C
243C     Call to the procedure for the calculations of the boundary conditions     
244      IF (present(procname)) THEN
245      Call Agrif_CorrectVariable
246     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname)
247      ELSE
248      Call Agrif_CorrectVariable
249     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight)
250      ENDIF
251C
252      child % var % oldvalues2D => childtemp % var % oldvalues2D
253      child % var % list_interp => childtemp % var %list_interp     
254C     
255      deallocate(childtemp % var)
256C
257C       
258      End Subroutine Agrif_Interp_bc_3D               
259C
260C
261C
262C
263C     **************************************************************************
264CCC   Subroutine Agrif_Interp_bc_4d
265C     **************************************************************************
266C 
267      Subroutine Agrif_Interp_bc_4d(TypeInterp,parent,child,tab,deb,fin,
268     &                             weight,pweight,procname)           
269C
270CCC   Description:
271CCC   Subroutine to calculate the boundary conditions on a fine grid for a 4D 
272CCC   grid variable.
273C
274C     Declarations:
275C     
276     
277C
278C     Arguments     
279      External :: procname
280      Optional :: procname     
281      INTEGER,DIMENSION(6) :: TypeInterp      ! TYPE of interpolation (linear,
282                                              ! lagrange, spline, ... )
283      TYPE(AGRIF_PVariable) :: parent         ! Variable on the parent grid
284      TYPE(AGRIF_PVariable) :: child          ! Variable on the child grid
285      TYPE(AGRIF_PVariable) :: childtemp      ! Temporary varaiable on the child
286                                              ! grid
287      INTEGER :: deb,fin                      ! Positions where interpolations
288                                              ! are done on the fine grid
289      REAL, DIMENSION(
290     &        lbound(child%var%array4,1):ubound(child%var%array4,1),
291     &        lbound(child%var%array4,2):ubound(child%var%array4,2),
292     &        lbound(child%var%array4,3):ubound(child%var%array4,3),
293     &        lbound(child%var%array4,4):ubound(child%var%array4,4)
294     &        ), Target :: tab ! Values of the grid variable
295      LOGICAL :: pweight                      ! Indicates if weight is used for
296                                              ! the temporal interpolation
297      REAL :: weight                          ! Coefficient for the time
298                                              ! interpolation
299C
300C
301C     Definition of a temporary AGRIF_PVariable data TYPE representing the grid 
302C     variable. 
303C
304      allocate(childtemp % var) 
305C
306      childtemp % var % root_var => child % var % root_var 
307C     
308C     Values of the grid variable
309      childtemp % var % array4 => tab 
310C
311C     Temporary results for the time interpolation before and after the space 
312C     interpolation       
313      childtemp % var % oldvalues2D => child % var % oldvalues2D
314C 
315C     Index indicating if a space interpolation is necessary
316      childtemp % var % interpIndex => child % var % interpIndex       
317      childtemp % var % Interpolationshouldbemade = 
318     &                 child % var % Interpolationshouldbemade 
319      childtemp % var % list_interp => child % var% list_interp         
320C
321C     Call to the procedure for the calculations of the boundary conditions
322      IF (present(procname)) THEN
323      Call Agrif_CorrectVariable
324     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname)
325      ELSE
326      Call Agrif_CorrectVariable
327     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight)
328      ENDIF
329C
330      child % var % oldvalues2D => childtemp % var % oldvalues2D
331      child % var % list_interp => childtemp % var %list_interp     
332C     
333      deallocate(childtemp % var)
334C
335C       
336      End Subroutine Agrif_Interp_bc_4D
337C
338C
339C
340C     **************************************************************************
341CCC   Subroutine Agrif_Interp_bc_5d
342C     **************************************************************************
343C 
344      Subroutine Agrif_Interp_bc_5d(TypeInterp,parent,child,tab,deb,fin,
345     &                             weight,pweight,procname)           
346C
347CCC   Description:
348CCC   Subroutine to calculate the boundary conditions on a fine grid for a 5D 
349CCC   grid variable.
350C
351C     Declarations:
352C     
353     
354C
355C     Arguments     
356      External :: procname
357      Optional :: procname     
358      INTEGER,DIMENSION(6) :: TypeInterp      ! TYPE of interpolation (linear,
359                                              ! lagrange, spline, ... )
360      TYPE(AGRIF_PVariable) :: parent         ! Variable on the parent grid
361      TYPE(AGRIF_PVariable) :: child          ! Variable on the child grid
362      TYPE(AGRIF_PVariable) :: childtemp      ! Temporary varaiable on the child
363                                              ! grid
364      INTEGER :: deb,fin                      ! Positions where interpolations
365                                              ! are done on the fine grid
366      REAL, DIMENSION(
367     &         lbound(child%var%array5,1):ubound(child%var%array5,1),
368     &         lbound(child%var%array5,2):ubound(child%var%array5,2),
369     &         lbound(child%var%array5,3):ubound(child%var%array5,3),
370     &         lbound(child%var%array5,4):ubound(child%var%array5,4),
371     &         lbound(child%var%array5,5):ubound(child%var%array5,5)
372     &         ), Target :: tab ! Values of the grid variable
373      LOGICAL :: pweight                      ! Indicates if weight is used for
374                                              ! the temporal interpolation
375      REAL :: weight                          ! Coefficient for the time
376                                              ! interpolation
377C
378C
379C     Definition of a temporary AGRIF_PVariable data TYPE representing the grid 
380C     variable. 
381C
382      allocate(childtemp % var) 
383C
384      childtemp % var % root_var => child % var % root_var 
385C     
386C     Values of the grid variable
387      childtemp % var % array5 => tab 
388C
389C     Temporary results for the time interpolation before and after the space 
390C     interpolation       
391      childtemp % var % oldvalues2D => child % var % oldvalues2D
392C 
393C     Index indicating if a space interpolation is necessary
394      childtemp % var % interpIndex => child % var % interpIndex       
395      childtemp % var % Interpolationshouldbemade = 
396     &                 child % var % Interpolationshouldbemade 
397      childtemp % var % list_interp => child % var% list_interp           
398C
399C     Call to the procedure for the calculations of the boundary conditions 
400      IF (present(procname)) THEN
401      Call Agrif_CorrectVariable
402     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname)
403      ELSE
404      Call Agrif_CorrectVariable
405     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight)
406      ENDIF
407C
408      child % var % oldvalues2D => childtemp % var % oldvalues2D
409      child % var % list_interp => childtemp % var %list_interp     
410C     
411      deallocate(childtemp % var)
412C
413C       
414      End Subroutine Agrif_Interp_bc_5D
415C
416C
417C
418C
419C     **************************************************************************
420CCC   Subroutine Agrif_Interp_bc_6d
421C     **************************************************************************
422C 
423      Subroutine Agrif_Interp_bc_6d(TypeInterp,parent,child,tab,deb,fin,
424     &                             weight,pweight)           
425C
426CCC   Description:
427CCC   Subroutine to calculate the boundary conditions on a fine grid for a 6D 
428CCC   grid variable.
429C
430C     Declarations:
431C     
432     
433C
434C     Arguments     
435      INTEGER,DIMENSION(6) :: TypeInterp      ! TYPE of interpolation (linear,
436                                              ! lagrange, spline, ... )
437      TYPE(AGRIF_PVariable) :: parent         ! Variable on the parent grid
438      TYPE(AGRIF_PVariable) :: child          ! Variable on the child grid
439      TYPE(AGRIF_PVariable) :: childtemp      ! Temporary varaiable on the child
440                                              ! grid
441      INTEGER :: deb,fin                      ! Positions where interpolations
442                                              ! are done on the fine grid
443      REAL, DIMENSION(
444     &         lbound(child%var%array6,1):ubound(child%var%array6,1),
445     &         lbound(child%var%array6,2):ubound(child%var%array6,2),
446     &         lbound(child%var%array6,3):ubound(child%var%array6,3),
447     &         lbound(child%var%array6,4):ubound(child%var%array6,4),
448     &         lbound(child%var%array6,5):ubound(child%var%array6,5),
449     &         lbound(child%var%array6,6):ubound(child%var%array6,6)
450     &         ), Target :: tab ! Values of the grid variable
451      LOGICAL :: pweight                      ! Indicates if weight is used for
452                                              ! the temporal interpolation
453      REAL :: weight                          ! Coefficient for the time
454                                              ! interpolation
455C
456C
457C     Definition of a temporary AGRIF_PVariable data TYPE representing the grid 
458C     variable. 
459C
460      allocate(childtemp % var) 
461C
462      childtemp % var % root_var => child % var % root_var 
463C     
464C     Values of the grid variable
465      childtemp % var % array6 => tab 
466C
467C     Temporary results for the time interpolation before and after the space 
468C     interpolation       
469      childtemp % var % oldvalues2D => child % var % oldvalues2D
470C 
471C     Index indicating if a space interpolation is necessary
472      childtemp % var % interpIndex => child % var % interpIndex       
473      childtemp % var % Interpolationshouldbemade = 
474     &                 child % var % Interpolationshouldbemade 
475      childtemp % var % list_interp => child % var% list_interp           
476C
477C     Call to the procedure for the calculations of the boundary conditions
478      Call Agrif_CorrectVariable
479     &     (TypeInterp,parent,childtemp,deb,fin,pweight,weight)
480C
481      child % var % oldvalues2D => childtemp % var % oldvalues2D
482      child % var % list_interp => childtemp % var %list_interp     
483C     
484      deallocate(childtemp % var)
485C
486C       
487      End Subroutine Agrif_Interp_bc_6D
488C
489C
490C     **************************************************************************
491CCC   Subroutine Agrif_CorrectVariable
492C     **************************************************************************
493C
494      Subroutine AGRIF_CorrectVariable(TypeInterp,parent,child,deb,fin,
495     &                                 pweight,weight,procname)
496C
497CCC   Description:
498CCC   Subroutine to calculate the boundary conditions on a fine grid.
499C
500C     Declarations:
501C     
502     
503C
504C     Arguments
505      External :: procname
506      Optional :: procname
507      TYPE(AGRIF_PVariable) :: parent         ! Variable on the parent grid
508      TYPE(AGRIF_PVariable) :: child          ! Variable on the child grid
509      INTEGER,DIMENSION(6)  :: TypeInterp     ! TYPE of interpolation
510                                              !    (linear,lagrange,...)
511      INTEGER               :: deb,fin        ! Positions where boundary
512                                              !    conditions are calculated
513      LOGICAL               :: pweight        ! Indicates if weight is used
514                                              !    for the time interpolation
515      REAL                  :: weight         ! Coefficient for the time
516                                              !    interpolation
517C
518C     Local scalars
519      TYPE(Agrif_Grid)    , Pointer :: Agrif_Child_Gr,Agrif_Parent_Gr
520      TYPE(AGRIF_Variable), Pointer :: root   ! Variable on the root grid
521      INTEGER                       :: nbdim  ! Number of dimensions of
522                                              !    the grid variable
523      INTEGER                       :: n
524      INTEGER,DIMENSION(6)          :: pttab_child  ! Index of the first point
525                                                    !    inside the domain for
526                                                    !    the child grid variable
527      INTEGER,DIMENSION(6)          :: pttab_parent ! Index of the first point
528                                                    !    inside the domain for
529                                                    !    the parent grid
530                                                    !    variable
531      INTEGER,DIMENSION(6)          :: nbtab_Child  ! Number of the cells
532      INTEGER,DIMENSION(6)          :: posvartab_Child    ! Position of the
533                                                    !    variable on the cell
534      INTEGER,DIMENSION(6)          :: loctab_Child ! Indicates if the child
535                                                    !    grid has a common
536                                                    !    border with the root
537                                                    !    grid     
538      REAL, DIMENSION(6)            :: s_child,s_parent   ! Positions of the
539                                                    !    parent and child grids
540      REAL, DIMENSION(6)            :: ds_child,ds_parent ! Space steps of the
541                                                    !    parent and child grids
542C
543C
544      loctab_child(:) = 0
545C
546      Agrif_Child_Gr => Agrif_Curgrid
547      Agrif_Parent_Gr => Agrif_Curgrid % parent
548      root => child % var % root_var 
549      nbdim = root % nbdim
550C
551      do n = 1,nbdim
552        posvartab_child(n) = root % posvar(n)     
553      enddo
554C
555C
556      do n = 1,nbdim
557C
558        Select case(root % interptab(n))
559C
560          case('x') ! x DIMENSION
561C
562            nbtab_Child(n) = Agrif_Child_Gr % nb(1)
563            pttab_Child(n) = root % point(1)
564            pttab_Parent(n) = root % point(1)
565            s_Child(n) = Agrif_Child_Gr % Agrif_x(1)
566            s_Parent(n) = Agrif_Parent_Gr % Agrif_x(1)
567            ds_Child(n) = Agrif_Child_Gr % Agrif_d(1)
568            ds_Parent(n) = Agrif_Parent_Gr % Agrif_d(1)
569            if (root % posvar(n) == 2) then
570                s_Child(n) = s_Child(n) + ds_Child(n)/2.
571                s_Parent(n) = s_Parent(n) + ds_Parent(n)/2.
572            endif
573C
574            if (Agrif_CURGRID % NearRootBorder(1)) 
575     &         loctab_child(n) = -1
576            if (Agrif_CURGRID % DistantRootBorder(1))
577     &         loctab_child(n) = -2
578            if ((Agrif_CURGRID % NearRootBorder(1)) .AND. 
579     &          (Agrif_CURGRID % DistantRootBorder(1))) 
580     &         loctab_child(n) = -3
581C
582          case('y') ! y DIMENSION     
583C
584            nbtab_Child(n) = Agrif_Child_Gr % nb(2)
585            pttab_Child(n) = root % point(2)
586            pttab_Parent(n) = root % point(2)
587            s_Child(n) = Agrif_Child_Gr % Agrif_x(2)
588            s_Parent(n) = Agrif_Parent_Gr % Agrif_x(2) 
589            ds_Child(n) = Agrif_Child_Gr % Agrif_d(2)
590            ds_Parent(n) = Agrif_Parent_Gr % Agrif_d(2)
591            if (root % posvar(n) == 2) then       
592                s_Child(n) = s_Child(n) + ds_Child(n)/2.
593                s_Parent(n) = s_Parent(n) + ds_Parent(n)/2.
594            endif       
595C
596            if (Agrif_CURGRID % NearRootBorder(2)) 
597     &         loctab_child(n) = -1
598            if (Agrif_CURGRID % DistantRootBorder(2)) 
599     &         loctab_child(n) = -2
600            if ((Agrif_CURGRID % NearRootBorder(2)) .AND. 
601     &          (Agrif_CURGRID % DistantRootBorder(2))) 
602     &         loctab_child(n) = -3
603C
604          case('z') ! z DIMENSION
605C
606            nbtab_Child(n) = Agrif_Child_Gr % nb(3)
607            pttab_Child(n) = root % point(3)
608            pttab_Parent(n) = root % point(3)
609            s_Child(n) = Agrif_Child_Gr % Agrif_x(3)
610            s_Parent(n) = Agrif_Parent_Gr % Agrif_x(3)
611            ds_Child(n) = Agrif_Child_Gr % Agrif_d(3)
612            ds_Parent(n) = Agrif_Parent_Gr % Agrif_d(3)
613            if (root % posvar(n) == 2) then       
614                s_Child(n) = s_Child(n) + ds_Child(n)/2.
615                s_Parent(n) = s_Parent(n) + ds_Parent(n)/2.
616            endif       
617C
618            if (Agrif_CURGRID % NearRootBorder(3)) 
619     &         loctab_child(n) = -1
620            if (Agrif_CURGRID % DistantRootBorder(3)) 
621     &         loctab_child(n) = -2
622            if ((Agrif_CURGRID % NearRootBorder(3)) .AND. 
623     &          (Agrif_CURGRID % DistantRootBorder(3))) 
624     &         loctab_child(n) = -3
625C
626          case('N') ! No space DIMENSION     
627C
628            select case (nbdim) 
629C     
630              case(1)
631                nbtab_Child(n) = SIZE(child % var % array1,n) - 1
632                pttab_Child(n) = lbound(child % var % array1,n)
633              case(2)
634                nbtab_Child(n) = SIZE(child % var % array2,n) - 1
635                pttab_Child(n) = lbound(child % var % array2,n)
636              case(3)
637                nbtab_Child(n) = SIZE(child % var % array3,n) - 1
638                pttab_Child(n) = lbound(child % var % array3,n) 
639              case(4)
640                nbtab_Child(n) = SIZE(child % var % array4,n) - 1
641                pttab_Child(n) = lbound(child % var % array4,n)
642              case(5)
643                nbtab_Child(n) = SIZE(child % var % array5,n) - 1
644                pttab_Child(n) = lbound(child % var % array5,n)     
645              case(6)
646                nbtab_Child(n) = SIZE(child % var % array6,n) - 1
647                pttab_Child(n) = lbound(child % var % array6,n)     
648C
649            end select
650C
651C           No interpolation but only a copy of the values of the grid variable
652C     
653            posvartab_child(n) = 1
654            pttab_Parent(n)= pttab_Child(n)
655            s_Child(n)=0.
656            s_Parent(n)=0. 
657            ds_Child(n)=1.
658            ds_Parent(n)=1.
659            loctab_child(n) = -3
660C
661        End select
662C
663      enddo
664C
665         IF (present(procname)) THEN
666          Call AGRIF_CorrectnD
667     &         (TypeInterp,parent,child,deb,fin,pweight,weight,
668     &          pttab_Child(1:nbdim),pttab_Parent(1:nbdim),
669     &          nbtab_Child(1:nbdim),posvartab_Child(1:nbdim),
670     &          loctab_Child(1:nbdim),
671     &          s_Child(1:nbdim),s_Parent(1:nbdim),
672     &          ds_Child(1:nbdim),ds_Parent(1:nbdim),nbdim,procname)
673        ELSE
674          Call AGRIF_CorrectnD
675     &         (TypeInterp,parent,child,deb,fin,pweight,weight,
676     &          pttab_Child(1:nbdim),pttab_Parent(1:nbdim),
677     &          nbtab_Child(1:nbdim),posvartab_Child(1:nbdim),
678     &          loctab_Child(1:nbdim),
679     &          s_Child(1:nbdim),s_Parent(1:nbdim),
680     &          ds_Child(1:nbdim),ds_Parent(1:nbdim),nbdim)
681         ENDIF
682C
683C
684      End subroutine AGRIF_CorrectVariable
685C
686C     **************************************************************************
687CCC   Subroutine Agrif_Correctnd
688C     **************************************************************************
689C
690      Subroutine AGRIF_Correctnd(TypeInterp,parent,child,deb,fin,
691     &                           pweight,weight,
692     &                           pttab_child,pttab_Parent,
693     &                           nbtab_Child,posvartab_Child,
694     &                           loctab_Child,
695     &                           s_Child,s_Parent,
696     &                           ds_Child,ds_Parent,nbdim,procname)
697C
698CCC   Description:
699CCC   Subroutine to calculate the boundary conditions for a nD grid variable on 
700CCC   a fine grid by using a space and time interpolations; it is called by the 
701CCC   Agrif_CorrectVariable procedure.
702C
703C
704C     Declarations:
705C
706     
707C
708#ifdef AGRIF_MPI
709C
710#include "mpif.h"
711C
712#endif
713C
714C     Arguments
715      External :: procname
716      Optional :: procname
717      INTEGER,DIMENSION(6) :: TypeInterp ! TYPE of interpolation (linear,
718                                         !   spline,...) 
719      TYPE(AGRIF_PVariable)    :: parent ! Variable on the parent grid
720      TYPE(AGRIF_PVariable)    :: child  ! Variable on the child grid
721      INTEGER                  :: deb,fin ! Positions where interpolations
722                                         !    are done
723      LOGICAL                  :: pweight ! Indicates if weight is used for
724                                         !    the temporal interpolation
725      REAL                     :: weight ! Coefficient for the temporal
726                                         !    interpolation
727      INTEGER                  :: nbdim  ! Number of dimensions of the grid
728                                         !    variable
729      INTEGER,DIMENSION(nbdim) :: pttab_child ! Index of the first point inside
730                                         !    the domain for the parent
731                                         !    grid variable
732      INTEGER,DIMENSION(nbdim) :: pttab_Parent ! Index of the first point
733                                         !   inside the domain for the
734                                         !   child grid variable
735      INTEGER,DIMENSION(nbdim) :: nbtab_Child ! Number of cells of the child
736                                         !    grid
737      INTEGER,DIMENSION(nbdim) :: posvartab_Child ! Position of the grid
738                                         !    variable (1 or 2)
739      INTEGER,DIMENSION(nbdim) :: loctab_Child ! Indicates if the child
740                                        !    grid has a common border with
741                                        !    the root grid
742      REAL   ,DIMENSION(nbdim) :: s_Child,s_Parent ! Positions of the parent
743                                        !   and child grids
744      REAL   ,DIMENSION(nbdim) :: ds_Child,ds_Parent ! Space steps of the
745                                        !    parent and child grids
746C
747C     Local variables
748      TYPE(AGRIF_PVariable)        :: restore ! Variable on the parent     
749      INTEGER,DIMENSION(nbdim,2)   :: lubglob
750      INTEGER                      :: i         
751      INTEGER                      :: kindex  ! Index used for safeguard
752                                       !    and time interpolation
753      INTEGER,DIMENSION(nbdim,2,2) :: indtab ! Arrays indicating the limits
754                                       !    of the child     
755      INTEGER,DIMENSION(nbdim,2,2) :: indtruetab ! grid variable where
756                                       !   boundary conditions are
757      INTEGER,DIMENSION(nbdim,2,2,nbdim)   :: ptres,ptres2 ! calculated
758      INTEGER                      :: nb,ndir,n,sizetab(1)
759      REAL, DIMENSION(:), Allocatable :: tab ! Array used for the interpolation
760      REAL    :: c1t,c2t               ! Coefficients for the time interpolation
761                                       !    (c2t=1-c1t)
762C
763#ifdef AGRIF_MPI
764C
765      INTEGER,DIMENSION(nbdim)   :: lower,upper
766      INTEGER,DIMENSION(nbdim)   :: ltab,utab
767      INTEGER,DIMENSION(nbdim)   :: lb,ub
768      INTEGER,DIMENSION(nbdim,2) :: iminmaxg
769      INTEGER                    :: code
770C
771#endif           
772C     
773C
774      indtab(1:nbdim,2,1) = pttab_child(1:nbdim) + nbtab_child(1:nbdim)
775     &          + deb
776      indtab(1:nbdim,2,2) = indtab(1:nbdim,2,1) + ( fin - deb )
777       
778      indtab(1:nbdim,1,1) = pttab_child(1:nbdim) - fin
779      indtab(1:nbdim,1,2) = pttab_child(1:nbdim) - deb
780                 
781      WHERE (posvartab_child(1:nbdim) == 2)
782        indtab(1:nbdim,1,1) = indtab(1:nbdim,1,1) - 1
783        indtab(1:nbdim,1,2) = indtab(1:nbdim,1,2) - 1
784      END WHERE
785     
786
787#if !defined AGRIF_MPI
788      Call Agrif_nbdim_Get_bound_dimension(child%var,lubglob(:,1),
789     &              lubglob(:,2),nbdim)
790C
791#else
792C     
793      Call Agrif_nbdim_Get_bound_dimension(child%var,lb,ub,nbdim)
794     
795      DO i = 1,nbdim
796C
797        Call Agrif_Invloc(lb(i),Agrif_Procrank,i,iminmaxg(i,1))
798        Call Agrif_Invloc(ub(i),Agrif_Procrank,i,iminmaxg(i,2))
799C
800      ENDDO
801C
802      iminmaxg(1:nbdim,2) = - iminmaxg(1:nbdim,2)
803     
804      CALL MPI_ALLREDUCE(iminmaxg,lubglob,2*nbdim,MPI_INTEGER,MPI_MIN,
805     &                     MPI_COMM_WORLD,code) 
806     
807      lubglob(1:nbdim,2) = - lubglob(1:nbdim,2)     
808C
809#endif
810C     
811      indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1),
812     &     lubglob(1:nbdim,1))
813      indtruetab(1:nbdim,1,2) = max(indtab(1:nbdim,1,2),
814     &     lubglob(1:nbdim,1))
815      indtruetab(1:nbdim,2,1) = min(indtab(1:nbdim,2,1),
816     &     lubglob(1:nbdim,2))
817      indtruetab(1:nbdim,2,2) = min(indtab(1:nbdim,2,2),
818     &     lubglob(1:nbdim,2))
819
820                             
821C 
822C
823      do nb = 1,nbdim
824C
825        do ndir = 1,2     
826C
827          if (loctab_child(nb) /= (-ndir) 
828     &        .AND. loctab_child(nb) /= -3) then
829C           
830              do n = 1,2
831C
832                ptres(nb,n,ndir,nb) = indtruetab(nb,ndir,n) 
833C
834              enddo             
835C
836              do i = 1,nbdim
837C     
838                if (i .NE. nb) then     
839C
840                    if (loctab_child(i) == -1 
841     &                                 .OR. loctab_child(i) == -3) then
842C
843                        ptres(i,1,ndir,nb) = pttab_child(i)
844C
845                      else
846C
847                        ptres(i,1,ndir,nb) = indtruetab(i,1,1)
848C
849                    endif
850C
851                    if (loctab_child(i) == -2 
852     &                                 .OR. loctab_child(i) == -3) then
853C
854                        if (posvartab_child(i) == 1) then
855C
856                            ptres(i,2,ndir,nb) = pttab_child(i)
857     &                                                 + nbtab_child(i)
858C
859                          else
860C
861                            ptres(i,2,ndir,nb) = pttab_child(i)
862     &                                             + nbtab_child(i) - 1
863C
864                        endif                             
865C
866                      else
867C
868                        ptres(i,2,ndir,nb) = indtruetab(i,2,2)
869C
870                    endif                       
871C     
872                endif
873C     
874              enddo
875     
876C
877#if defined AGRIF_MPI
878                 Call Agrif_nbdim_Get_bound_dimension
879     &                    (child%var,lower,upper,nbdim)
880
881              do i = 1,nbdim
882C
883                 Call GetLocalBoundaries(ptres(i,1,ndir,nb),
884     &                                  ptres(i,2,ndir,nb),i,
885     &                                  lower(i),upper(i),
886     &                                  ltab(i),utab(i))
887                 ptres2(i,1,ndir,nb) = max(ltab(i),lower(i)) 
888                 ptres2(i,2,ndir,nb) = min(utab(i),upper(i))
889                 if ((i == nb) .AND. (ndir == 1)) then
890                     ptres2(i,2,ndir,nb) = max(utab(i),lower(i))
891                 elseif ((i == nb) .AND. (ndir == 2)) then
892                     ptres2(i,1,ndir,nb) = min(ltab(i),upper(i)) 
893                 endif
894C
895              enddo
896#else
897              ptres2(:,:,ndir,nb) = ptres(:,:,ndir,nb)
898#endif     
899           
900        endif
901     
902        enddo
903       enddo
904C
905      if (child % var % interpIndex 
906     &        /= Agrif_Curgrid % parent % ngridstep .OR.
907     &    child%var%Interpolationshouldbemade ) then
908C
909C     Space interpolation 
910C     
911      kindex = 1   
912C
913      do nb = 1,nbdim
914C
915        do ndir = 1,2               
916C
917          if (loctab_child(nb) /= (-ndir) 
918     &        .AND. loctab_child(nb) /= -3) then 
919C
920              IF (present(procname)) THEN
921              Call Agrif_InterpnD
922     &             (TYPEInterp,parent,child,
923     &              ptres(1:nbdim,1,ndir,nb),ptres(1:nbdim,2,ndir,nb),
924     &              pttab_child(1:nbdim),pttab_Parent(1:nbdim),
925     &              s_Child(1:nbdim),s_Parent(1:nbdim),
926     &              ds_Child(1:nbdim),ds_Parent(1:nbdim),
927     &              restore,.FALSE.,nbdim,procname)
928              ELSE
929              Call Agrif_InterpnD             
930     &             (TYPEInterp,parent,child,
931     &              ptres(1:nbdim,1,ndir,nb),ptres(1:nbdim,2,ndir,nb),
932     &              pttab_child(1:nbdim),pttab_Parent(1:nbdim),
933     &              s_Child(1:nbdim),s_Parent(1:nbdim),
934     &              ds_Child(1:nbdim),ds_Parent(1:nbdim),
935     &              restore,.FALSE.,nbdim)
936              ENDIF
937     
938              IF (.NOT. child%var%interpolationshouldbemade) THEN
939C     
940C             Safeguard of the values of the grid variable (at times n and n+1
941C                on the parent grid)
942C     
943              sizetab(1) = 1
944C
945              do i = 1,nbdim
946C         
947                sizetab(1) = sizetab(1)
948     &              * (ptres2(i,2,ndir,nb)-ptres2(i,1,ndir,nb)+1)
949C     
950              enddo
951
952              Call saveAfterInterp(child,
953     &           ptres2(:,:,ndir,nb),kindex,sizetab(1),nbdim)
954C
955           ENDIF
956C
957          endif
958C
959        enddo       
960C
961      enddo     
962C
963C
964      child % var % interpIndex = Agrif_Curgrid % parent % ngridstep
965C
966C
967      endif 
968     
969              IF (.NOT. child%var%interpolationshouldbemade) THEN
970C
971C
972C     Calculation of the coefficients c1t and c2t for the temporary
973C        interpolation
974      if (pweight) then
975C
976          c1t = weight
977C
978        else
979C
980          c1t = (REAL(AGRIF_Nbstepint()) + 1.) / Agrif_Rhot()
981C
982      endif
983C                                 
984      c2t = 1. - c1t
985C           
986C     Time interpolation
987C
988      kindex = 1
989C
990      do nb = 1,nbdim
991C
992        do ndir = 1,2     
993C
994          if (loctab_child(nb) /= (-ndir) 
995     &        .AND. loctab_child(nb) /= -3) then
996
997              Call timeInterpolation
998     &             (child,ptres2(:,:,ndir,nb),kindex,c1t,c2t,nbdim)
999          endif
1000C
1001        enddo
1002C     
1003      enddo
1004C     
1005
1006       ENDIF
1007C 
1008      End Subroutine Agrif_Correctnd
1009C
1010C
1011C     **************************************************************************
1012CCC   Subroutine saveAfterInterp
1013C     **************************************************************************
1014C
1015      Subroutine saveAfterInterp(child,bounds,kindex,newsize,nbdim)
1016C
1017CCC   Descritpion:
1018CCC   Subroutine used to save the values of the grid variable on the fine grid 
1019CCC   after the space interpolation. 
1020C
1021C     Declarations:
1022C
1023     
1024C
1025C     argument
1026      TYPE (AGRIF_PVariable) :: child   ! The fine grid variable
1027      INTEGER                :: kindex  ! Index indicating where this safeguard
1028                                        ! is done on the fine grid
1029      INTEGER :: nbdim, newsize
1030      INTEGER,DIMENSION(nbdim,2) :: bounds
1031C
1032C     Local scalars
1033      INTEGER                      :: ir,jr,kr,lr,mr,nr
1034C
1035C
1036C     Allocation of the array oldvalues2d
1037
1038C
1039      if (newsize .LE. 0) return
1040C
1041      Call Agrif_Checksize
1042     &     (child,kindex+newsize) 
1043
1044      if (child % var % interpIndex 
1045     &        /= Agrif_Curgrid % parent % ngridstep ) then
1046       child%var%oldvalues2d(kindex:kindex+newsize-1,1)=
1047     &        child%var%oldvalues2d(kindex:kindex+newsize-1,2)
1048      endif
1049
1050      SELECT CASE (nbdim)
1051      CASE (1)
1052
1053         do ir=bounds(1,1),bounds(1,2)
1054            child%var%oldvalues2d(kindex,2) =
1055     &           child%var%array1(ir)
1056            kindex = kindex + 1
1057         enddo       
1058C
1059      CASE (2)
1060
1061        do jr=bounds(2,1),bounds(2,2)
1062           do ir=bounds(1,1),bounds(1,2)
1063            child%var%oldvalues2d(kindex,2) =
1064     &           child%var%array2(ir,jr)
1065            kindex = kindex + 1
1066           enddo
1067         enddo     
1068C
1069      CASE (3)
1070        do kr=bounds(3,1),bounds(3,2)
1071          do jr=bounds(2,1),bounds(2,2)
1072             do ir=bounds(1,1),bounds(1,2)
1073            child%var%oldvalues2d(kindex,2) =
1074     &           child%var%array3(ir,jr,kr)
1075            kindex = kindex + 1
1076             enddo
1077           enddo
1078         enddo       
1079C
1080      CASE (4)
1081          do lr=bounds(4,1),bounds(4,2)
1082           do kr=bounds(3,1),bounds(3,2)
1083             do jr=bounds(2,1),bounds(2,2)
1084               do ir=bounds(1,1),bounds(1,2)
1085            child%var%oldvalues2d(kindex,2) =
1086     &           child%var%array4(ir,jr,kr,lr)
1087            kindex = kindex + 1
1088               enddo
1089             enddo
1090           enddo
1091         enddo         
1092C
1093      CASE (5)
1094         do mr=bounds(5,1),bounds(5,2)
1095           do lr=bounds(4,1),bounds(4,2)
1096            do kr=bounds(3,1),bounds(3,2) 
1097              do jr=bounds(2,1),bounds(2,2)
1098                 do ir=bounds(1,1),bounds(1,2)
1099            child%var%oldvalues2d(kindex,2) =
1100     &           child%var%array5(ir,jr,kr,lr,mr)
1101            kindex = kindex + 1
1102                 enddo
1103               enddo
1104             enddo
1105           enddo
1106         enddo   
1107C
1108      CASE (6)
1109        do nr=bounds(6,1),bounds(6,2)
1110           do mr=bounds(5,1),bounds(5,2)
1111             do lr=bounds(4,1),bounds(4,2)
1112               do kr=bounds(3,1),bounds(3,2)
1113                do jr=bounds(2,1),bounds(2,2)
1114                   do ir=bounds(1,1),bounds(1,2)
1115            child%var%oldvalues2d(kindex,2) =
1116     &           child%var%array6(ir,jr,kr,lr,mr,nr)
1117            kindex = kindex + 1
1118                   enddo
1119                 enddo
1120               enddo
1121             enddo
1122           enddo
1123         enddo       
1124      END SELECT
1125C
1126C                                                 
1127      End subroutine saveAfterInterp
1128C
1129C
1130C
1131C     **************************************************************************
1132CCC   Subroutine timeInterpolation
1133C     **************************************************************************
1134C
1135      Subroutine timeInterpolation(child,bounds,kindex,c1t,c2t,nbdim) 
1136C
1137CCC   Descritpion:
1138CCC   Subroutine for a linear time interpolation on the child grid. 
1139C
1140C     Declarations:
1141C
1142     
1143C
1144C     argument
1145      TYPE (AGRIF_PVariable) :: child  ! The fine grid variable
1146      INTEGER :: nbdim
1147      INTEGER,DIMENSION(nbdim,2) :: bounds
1148      INTEGER                :: kindex ! Index indicating the values of the fine
1149                                       ! grid got before and after the space
1150                                       ! interpolation and used for the time
1151                                       ! interpolation
1152      REAL                   :: c1t,c2t! coefficients for the time interpolation
1153                                       ! (c2t=1-c1t) 
1154C
1155C     Local aruments     
1156      INTEGER :: i
1157C     Local scalars
1158      INTEGER                      :: ir,jr,kr,lr,mr,nr     
1159C
1160C
1161     
1162      SELECT CASE (nbdim)
1163      CASE (1)
1164
1165         do ir=bounds(1,1),bounds(1,2)
1166                child%var%array1(ir) =
1167     &           c2t*child % var % oldvalues2d(kindex,1)   
1168     &         + c1t*child % var % oldvalues2d(kindex,2)   
1169            kindex = kindex + 1
1170         enddo       
1171C
1172      CASE (2)
1173
1174        do jr=bounds(2,1),bounds(2,2)
1175           do ir=bounds(1,1),bounds(1,2)
1176                child%var%array2(ir,jr) =
1177     &           c2t*child % var % oldvalues2d(kindex,1)   
1178     &         + c1t*child % var % oldvalues2d(kindex,2) 
1179            kindex = kindex + 1
1180           enddo
1181         enddo     
1182C
1183      CASE (3)
1184        do kr=bounds(3,1),bounds(3,2)
1185          do jr=bounds(2,1),bounds(2,2)
1186             do ir=bounds(1,1),bounds(1,2)
1187                child%var%array3(ir,jr,kr) =
1188     &           c2t*child % var % oldvalues2d(kindex,1)   
1189     &         + c1t*child % var % oldvalues2d(kindex,2) 
1190            kindex = kindex + 1
1191             enddo
1192           enddo
1193         enddo       
1194C
1195      CASE (4)
1196          do lr=bounds(4,1),bounds(4,2)
1197           do kr=bounds(3,1),bounds(3,2)
1198             do jr=bounds(2,1),bounds(2,2)
1199               do ir=bounds(1,1),bounds(1,2)
1200                child%var%array4(ir,jr,kr,lr) =
1201     &           c2t*child % var % oldvalues2d(kindex,1)   
1202     &         + c1t*child % var % oldvalues2d(kindex,2) 
1203            kindex = kindex + 1
1204               enddo
1205             enddo
1206           enddo
1207         enddo         
1208C
1209      CASE (5)
1210         do mr=bounds(5,1),bounds(5,2)
1211           do lr=bounds(4,1),bounds(4,2)
1212            do kr=bounds(3,1),bounds(3,2) 
1213              do jr=bounds(2,1),bounds(2,2)
1214                 do ir=bounds(1,1),bounds(1,2)
1215                child%var%array5(ir,jr,kr,lr,mr) =
1216     &           c2t*child % var % oldvalues2d(kindex,1)   
1217     &         + c1t*child % var % oldvalues2d(kindex,2) 
1218            kindex = kindex + 1
1219                 enddo
1220               enddo
1221             enddo
1222           enddo
1223         enddo   
1224C
1225      CASE (6)
1226        do nr=bounds(6,1),bounds(6,2)
1227           do mr=bounds(5,1),bounds(5,2)
1228             do lr=bounds(4,1),bounds(4,2)
1229               do kr=bounds(3,1),bounds(3,2)
1230                do jr=bounds(2,1),bounds(2,2)
1231                   do ir=bounds(1,1),bounds(1,2)
1232                child%var%array6(ir,jr,kr,lr,mr,nr) =
1233     &           c2t*child % var % oldvalues2d(kindex,1)   
1234     &         + c1t*child % var % oldvalues2d(kindex,2) 
1235            kindex = kindex + 1
1236                   enddo
1237                 enddo
1238               enddo
1239             enddo
1240           enddo
1241         enddo       
1242      END SELECT
1243                                                   
1244C
1245C
1246      End subroutine timeInterpolation     
1247C
1248C
1249C
1250C     **************************************************************************
1251CCC   Subroutine Agrif_Checksize
1252C     **************************************************************************
1253C
1254      Subroutine Agrif_Checksize(child,newsize)
1255C
1256CCC   Descritpion:
1257CCC   Subroutine used in the saveAfterInterp procedure to allocate the 
1258CCC   oldvalues2d array. 
1259C
1260C     Declarations:
1261C
1262     
1263C     
1264C     TYPE argument
1265      TYPE (AGRIF_PVariable) :: child  ! The fine grid variable
1266C 
1267C     Scalar arguments
1268      INTEGER :: newsize               ! Size of the domains where the boundary
1269                                       ! conditions are calculated 
1270C
1271C     Local arrays
1272      REAL, DIMENSION(:,:), Allocatable :: tempoldvalues ! Temporary array
1273C               
1274C
1275      if (.NOT. associated(child % var % oldvalues2d)) then
1276C
1277          allocate(child % var % oldvalues2d(newsize,2))
1278C 
1279          child % var % oldvalues2d=0.
1280C
1281        else
1282C
1283          if (SIZE(child % var % oldvalues2d,1) < newsize) then   
1284C
1285              allocate(tempoldvalues(SIZE(child % var %
1286     &                                    oldvalues2d,1),2))
1287C
1288              tempoldvalues = child % var % oldvalues2d
1289C
1290              deallocate(child % var % oldvalues2d)
1291C
1292              allocate(child % var % oldvalues2d(newsize,2))
1293C           
1294              child%var%oldvalues2d=0.
1295C
1296              child % var % oldvalues2d(1:SIZE(tempoldvalues,1),:) = 
1297     &        tempoldvalues(:,:)
1298C
1299              deallocate(tempoldvalues)
1300C
1301          endif
1302C
1303      endif                                                   
1304C
1305C
1306      End  Subroutine Agrif_Checksize
1307C
1308C
1309C
1310C     
1311      End Module AGRIF_boundary
1312
Note: See TracBrowser for help on using the repository browser.