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

Last change on this file since 396 was 396, checked in by opalod, 18 years ago

Initial revision

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