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.F90 in vendors/AGRIF/CMEMS_2020/AGRIF_FILES – NEMO

source: vendors/AGRIF/CMEMS_2020/AGRIF_FILES/modbc.F90 @ 10087

Last change on this file since 10087 was 10087, checked in by rblod, 6 years ago

update AGRIF library

  • Property svn:keywords set to Id
File size: 31.4 KB
Line 
1!
2! $Id$
3!
4!     AGRIF (Adaptive Grid Refinement In Fortran)
5!
6!     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7!                        Christophe Vouland (Christophe.Vouland@imag.fr)
8!
9!     This program is free software; you can redistribute it and/or modify
10!     it under the terms of the GNU General Public License as published by
11!     the Free Software Foundation; either version 2 of the License, or
12!     (at your option) any later version.
13!
14!     This program is distributed in the hope that it will be useful,
15!     but WITHOUT ANY WARRANTY; without even the implied warranty of
16!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17!     GNU General Public License for more details.
18!
19!     You should have received a copy of the GNU General Public License
20!     along with this program; if not, write to the Free Software
21!     Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
22!
23!
24!> Module Agrif_Boundary.
25!>
26!> Contains subroutines to calculate the boundary conditions on the child grids from their
27!> parent grids.
28!
29module Agrif_Boundary
30!
31    use Agrif_Interpolation
32!
33    implicit none
34    REAL,DIMENSION(:),ALLOCATABLE :: parray_temp
35!
36contains
37!
38!===================================================================================================
39!  subroutine Agrif_CorrectVariable
40!
41!> subroutine to calculate the boundary conditions on a fine grid
42!---------------------------------------------------------------------------------------------------
43subroutine Agrif_CorrectVariable ( parent, child, pweight, weight, procname )
44!---------------------------------------------------------------------------------------------------
45    type(Agrif_Variable), pointer       :: parent       !< Variable on the parent grid
46    type(Agrif_Variable), pointer       :: child        !< Variable on the child grid
47    logical                             :: pweight      !< Indicates if weight is used for the time interpolation
48    real                                :: weight       !< Coefficient for the time interpolation
49    procedure()                         :: procname     !< Data recovery procedure
50!
51    type(Agrif_Grid)    , pointer :: Agrif_Child_Gr, Agrif_Parent_Gr
52    type(Agrif_Variable), pointer :: root_var   ! Variable on the root grid
53    integer                :: nbdim  ! Number of dimensions of the grid variable
54    integer                :: n
55    integer, dimension(6)  :: lb_child       ! Index of the first point inside the domain for
56                                                !    the child grid variable
57    integer, dimension(6)  :: lb_parent      ! Index of the first point inside the domain for
58                                                !    the parent grid variable
59    integer, dimension(6)  :: ub_child     !  Upper bound on the child grid
60    integer, dimension(6)  :: nb_child       ! Number of cells for child
61    integer, dimension(6)  :: posvartab_child   ! Position of the variable on the cell
62    integer, dimension(6)  :: loctab_child      ! Indicates if the child grid has a common border
63                                                !    with the root grid
64    real(kind=8), dimension(6)     :: s_child, s_parent   ! Positions of the parent and child grids
65    real(kind=8), dimension(6)     :: ds_child, ds_parent ! Space steps of the parent and child grids
66!
67    call PreProcessToInterpOrUpdate( parent,   child,       &
68                                     nb_child, ub_child,    &
69                                     lb_child, lb_parent,   &
70                                      s_child,  s_parent,   &
71                                     ds_child, ds_parent, nbdim, interp=.true.)
72    root_var => child % root_var
73    Agrif_Child_Gr => Agrif_Curgrid
74    Agrif_Parent_Gr => Agrif_Curgrid % parent
75!
76    loctab_child(:) = 0
77    posvartab_child(1:nbdim) = root_var % posvar(1:nbdim)
78!
79    do n = 1,nbdim
80!
81        select case(root_var % interptab(n))
82!
83        case('x') ! x DIMENSION
84!
85            if (Agrif_Curgrid % NearRootBorder(1))      loctab_child(n) = -1
86            if (Agrif_Curgrid % DistantRootBorder(1))   loctab_child(n) = -2
87            if ((Agrif_Curgrid % NearRootBorder(1)) .AND. &
88                (Agrif_Curgrid % DistantRootBorder(1))) loctab_child(n) = -3
89!
90        case('y') ! y DIMENSION
91!
92            if (Agrif_Curgrid % NearRootBorder(2))      loctab_child(n) = -1
93            if (Agrif_Curgrid % DistantRootBorder(2))   loctab_child(n) = -2
94            if ((Agrif_Curgrid % NearRootBorder(2)) .AND. &
95                (Agrif_Curgrid % DistantRootBorder(2))) loctab_child(n) = -3
96!
97        case('z') ! z DIMENSION
98!
99            if (Agrif_Curgrid % NearRootBorder(3))      loctab_child(n) = -1
100            if (Agrif_Curgrid % DistantRootBorder(3))   loctab_child(n) = -2
101            if ((Agrif_Curgrid % NearRootBorder(3)) .AND. &
102                (Agrif_Curgrid % DistantRootBorder(3))) loctab_child(n) = -3
103!
104        case('N') ! No space DIMENSION
105!
106            posvartab_child(n) = 1
107            loctab_child(n) = -3
108!
109        end select
110!
111    enddo
112!
113    call Agrif_Correctnd(parent, child, pweight, weight,                    &
114                         lb_child(1:nbdim), lb_parent(1:nbdim),             &
115                         nb_child(1:nbdim), posvartab_child(1:nbdim),       &
116                         loctab_child(1:nbdim),                             &
117                         s_child(1:nbdim), s_parent(1:nbdim),               &
118                         ds_child(1:nbdim),ds_parent(1:nbdim), nbdim, procname )
119!---------------------------------------------------------------------------------------------------
120end subroutine Agrif_CorrectVariable
121!===================================================================================================
122!
123!===================================================================================================
124!  subroutine Agrif_Correctnd
125!
126!> calculates the boundary conditions for a nD grid variable on a fine grid by using
127!> a space and time interpolations; it is called by the #Agrif_CorrectVariable procedure
128!---------------------------------------------------------------------------------------------------
129subroutine Agrif_Correctnd ( parent, child, pweight, weight,                        &
130                             pttab_child, pttab_Parent,                             &
131                             nbtab_Child, posvartab_Child, loctab_Child,            &
132                             s_Child, s_Parent, ds_Child, ds_Parent,                &
133                             nbdim, procname )
134!---------------------------------------------------------------------------------------------------
135#if defined AGRIF_MPI
136    include 'mpif.h'
137#endif
138!
139    TYPE(Agrif_Variable), pointer       :: parent       !< Variable on the parent grid
140    TYPE(Agrif_Variable), pointer       :: child        !< Variable on the child grid
141    LOGICAL                             :: pweight      !< Indicates if weight is used for the temporal interpolation
142    REAL                                :: weight       !< Coefficient for the temporal interpolation
143    INTEGER, DIMENSION(nbdim)   :: pttab_child          !< Index of the first point inside the domain for the parent grid variable
144    INTEGER, DIMENSION(nbdim)   :: pttab_Parent         !< Index of the first point inside the domain for the child  grid variable
145    INTEGER, DIMENSION(nbdim)   :: nbtab_Child          !< Number of cells of the child grid
146    INTEGER, DIMENSION(nbdim)   :: posvartab_Child      !< Position of the grid variable (1 or 2)
147    INTEGER, DIMENSION(nbdim)   :: loctab_Child         !< Indicates if the child grid has a common border with the root grid
148    REAL(kind=8)   , DIMENSION(nbdim)   :: s_Child,  s_Parent   !< Positions of the parent and child grids
149    REAL(kind=8)   , DIMENSION(nbdim)   :: ds_Child, ds_Parent  !< Space steps of the parent and child grids
150    INTEGER                             :: nbdim        !< Number of dimensions of the grid variable
151    procedure()                         :: procname     !< Data recovery procedure
152!
153    INTEGER,DIMENSION(6)                :: type_interp     ! Type of interpolation (linear, spline,...)
154    INTEGER,DIMENSION(6,6)              :: type_interp_bc  ! Type of interpolation (linear, spline,...)
155    INTEGER,DIMENSION(nbdim,2,2)        :: childarray
156    INTEGER,DIMENSION(nbdim,2)          :: lubglob
157    INTEGER                             :: kindex       ! Index used for safeguard and time interpolation
158    INTEGER,DIMENSION(nbdim,2,2)        :: indtab       ! Arrays indicating the limits of the child
159    INTEGER,DIMENSION(nbdim,2,2)        :: indtruetab   ! grid variable where boundary conditions are
160    INTEGER,DIMENSION(nbdim,2,2,nbdim)  :: ptres,ptres2 ! calculated
161    INTEGER,DIMENSION(nbdim)            :: coords
162    INTEGER                             :: i, nb, ndir,j,k,l
163    INTEGER                             :: n, sizetab
164    INTEGER                             :: ibeg, iend
165    INTEGER                             :: i1,i2,j1,j2,k1,k2,l1,l2,m1,m2,n1,n2
166    REAL                                :: c1t,c2t      ! Coefficients for the time interpolation (c2t=1-c1t)
167    INTEGER :: isize
168    INTEGER :: kindex_2d(2,nbdim)
169
170#if defined AGRIF_MPI
171!
172    INTEGER, DIMENSION(nbdim)   :: lower, upper
173    INTEGER, DIMENSION(nbdim)   :: ltab, utab
174!
175#endif
176!
177    type_interp_bc = child % root_var % type_interp_bc
178    coords         = child % root_var % coords
179!
180    ibeg = child % bcinf
181    iend = child % bcsup
182!
183    indtab(1:nbdim,2,1) = pttab_child(1:nbdim) + nbtab_child(1:nbdim) + ibeg
184    indtab(1:nbdim,2,2) = indtab(1:nbdim,2,1) + ( iend - ibeg )
185
186    indtab(1:nbdim,1,1) = pttab_child(1:nbdim) - iend
187    indtab(1:nbdim,1,2) = pttab_child(1:nbdim) - ibeg
188
189    WHERE (posvartab_child(1:nbdim) == 2)
190        indtab(1:nbdim,1,1) = indtab(1:nbdim,1,1) - 1
191        indtab(1:nbdim,1,2) = indtab(1:nbdim,1,2) - 1
192    END WHERE
193!
194!   call Agrif_get_var_global_bounds(child,lubglob,nbdim)
195    lubglob = child%lubglob(1:nbdim,:)
196!
197    indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1), lubglob(1:nbdim,1))
198    indtruetab(1:nbdim,1,2) = max(indtab(1:nbdim,1,2), lubglob(1:nbdim,1))
199    indtruetab(1:nbdim,2,1) = min(indtab(1:nbdim,2,1), lubglob(1:nbdim,2))
200    indtruetab(1:nbdim,2,2) = min(indtab(1:nbdim,2,2), lubglob(1:nbdim,2))
201   
202!
203    do nb = 1,nbdim
204        do ndir = 1,2
205!
206            if (loctab_child(nb) /= (-ndir) .AND. loctab_child(nb) /= -3) then
207!
208                do n = 1,2
209                    ptres(nb,n,ndir,nb) = indtruetab(nb,ndir,n)
210                enddo
211!
212                do i = 1,nbdim
213!
214                    if (i /= nb) then
215!
216                        if (loctab_child(i) == -1 .OR. loctab_child(i) == -3) then
217                            ptres(i,1,ndir,nb) = pttab_child(i)
218                        else
219                            ptres(i,1,ndir,nb) = indtruetab(i,1,1)
220                        endif
221                        if (loctab_child(i) == -2 .OR. loctab_child(i) == -3) then
222                            if (posvartab_child(i) == 1) then
223                                ptres(i,2,ndir,nb) = pttab_child(i) + nbtab_child(i)
224                            else
225                                ptres(i,2,ndir,nb) = pttab_child(i) + nbtab_child(i) - 1
226                            endif
227                        else
228                            ptres(i,2,ndir,nb) = indtruetab(i,2,2)
229                        endif
230!
231                    endif
232!
233                enddo
234
235!
236#if defined AGRIF_MPI
237                call Agrif_get_var_bounds_array(child,lower,upper,nbdim)
238
239                do i = 1,nbdim
240!
241                    Call Agrif_GetLocalBoundaries(ptres(i,1,ndir,nb), ptres(i,2,ndir,nb),  &
242                                                  coords(i), lower(i), upper(i), ltab(i), utab(i) )
243                    ptres2(i,1,ndir,nb) = max(ltab(i),lower(i))
244                    ptres2(i,2,ndir,nb) = min(utab(i),upper(i))
245                    if ((i == nb) .AND. (ndir == 1)) then
246                        ptres2(i,2,ndir,nb) = max(utab(i),lower(i))
247                    elseif ((i == nb) .AND. (ndir == 2)) then
248                        ptres2(i,1,ndir,nb) = min(ltab(i),upper(i))
249                    endif
250!
251                enddo
252#else
253                ptres2(:,:,ndir,nb) = ptres(:,:,ndir,nb)
254#endif
255            endif
256!
257        enddo   ! ndir = 1,2
258    enddo       ! nb = 1,nbdim
259!
260    if ( child % interpIndex /= Agrif_Curgrid % parent % ngridstep .OR. &
261         child % Interpolationshouldbemade ) then
262!
263!     Space interpolation
264!
265        kindex = 1
266!
267        do nb = 1,nbdim
268
269            type_interp = type_interp_bc(nb,:)
270
271            do ndir = 1,2
272!
273                if (loctab_child(nb) /= (-ndir) .AND. loctab_child(nb) /= -3) then
274!
275
276                    call Agrif_InterpnD(type_interp, parent, child,             &
277                                        ptres(1:nbdim,1,ndir,nb),               &
278                                        ptres(1:nbdim,2,ndir,nb),               &
279                                        pttab_child(1:nbdim),                   &
280                                        pttab_Parent(1:nbdim),                  &
281                                        s_Child(1:nbdim), s_Parent(1:nbdim),    &
282                                        ds_Child(1:nbdim),ds_Parent(1:nbdim),   &
283                                        NULL(), .FALSE., nbdim,                 &
284                                        childarray,                             &
285                                        child%memberin(nb,ndir), .TRUE., procname, coords(nb), ndir)
286
287                    child % childarray(1:nbdim,:,:,nb,ndir) = childarray
288
289                    if (.not. child%interpolationshouldbemade) then
290!
291!                       Safeguard of the values of the grid variable (at times n and n+1 on the parent grid)
292!
293                        sizetab = 1
294                        do i = 1,nbdim
295                            sizetab = sizetab * (ptres2(i,2,ndir,nb)-ptres2(i,1,ndir,nb)+1)
296                        enddo
297
298                        call saveAfterInterp(child,ptres2(:,:,ndir,nb),kindex,sizetab,nbdim)
299!
300                    endif
301!
302                endif
303!
304            enddo   ! ndir = 1,2
305        enddo       ! nb = 1,nbdim
306!
307        child % interpIndex = Agrif_Curgrid % parent % ngridstep
308!
309    endif
310!
311    if (.not. child%interpolationshouldbemade) then
312!
313!       Calculation of the coefficients c1t and c2t for the temporary interpolation
314!
315        if (pweight) then
316            c1t = weight
317        else
318            c1t = (REAL(Agrif_Nbstepint()) + 1.) / Agrif_Rhot()
319        endif
320        c2t = 1. - c1t
321!
322!       Time interpolation
323!
324        kindex = 1
325!
326        do nb = 1,nbdim
327            do ndir = 1,2
328                kindex_2d(ndir,nb) = kindex
329                if ( (loctab_child(nb) /= (-ndir)) .AND. (loctab_child(nb) /= -3) .AND. child%memberin(nb,ndir) ) then
330                    Call timeInterpolation(child,ptres2(:,:,ndir,nb),kindex,c1t,c2t,nbdim)
331                endif
332            enddo
333        enddo
334!
335    do nb = 1,nbdim
336    do ndir = 1,2
337        if ( (loctab_child(nb) /= (-ndir)) .AND. (loctab_child(nb) /= -3) .AND. child%memberin(nb,ndir) ) then
338
339         do i=1,nbdim
340         if (ptres2(i,1,ndir,nb) /= child%childarray(i,1,2,nb,ndir)) then
341            print *,'problem ptres2 childarray 1 ',ptres2(i,1,ndir,nb) /= child%childarray(i,1,2,nb,ndir)
342            stop
343         endif
344         if (ptres2(i,2,ndir,nb) /= child%childarray(i,2,2,nb,ndir)) then
345            print *,'problem ptres2 childarray 2 ',ptres2(i,2,ndir,nb) /= child%childarray(i,2,2,nb,ndir)
346          stop
347         endif
348         enddo
349
350            select case(nbdim)
351            case(1)
352                i1 = child % childarray(1,1,2,nb,ndir)
353                i2 = child % childarray(1,2,2,nb,ndir)
354
355                call procname(parray1(i1:i2),                               &
356                              i1,i2, .FALSE.,coords(nb),ndir)
357            case(2)
358                i1 = child % childarray(1,1,2,nb,ndir)
359                i2 = child % childarray(1,2,2,nb,ndir)
360                j1 = child % childarray(2,1,2,nb,ndir)
361                j2 = child % childarray(2,2,2,nb,ndir)
362
363                call procname(parray2(i1:i2,j1:j2),                         &
364                              i1,i2,j1,j2, .FALSE.,coords(nb),ndir)
365            case(3)
366
367                i1 = child % childarray(1,1,2,nb,ndir)
368                i2 = child % childarray(1,2,2,nb,ndir)
369                j1 = child % childarray(2,1,2,nb,ndir)
370                j2 = child % childarray(2,2,2,nb,ndir)
371                k1 = child % childarray(3,1,2,nb,ndir)
372                k2 = child % childarray(3,2,2,nb,ndir)
373
374               call procname(parray_temp(kindex_2d(ndir,nb)),i1,i2,j1,j2,k1,k2, .FALSE.,coords(nb),ndir)
375
376            case(4)
377                i1 = child % childarray(1,1,2,nb,ndir)
378                i2 = child % childarray(1,2,2,nb,ndir)
379                j1 = child % childarray(2,1,2,nb,ndir)
380                j2 = child % childarray(2,2,2,nb,ndir)
381                k1 = child % childarray(3,1,2,nb,ndir)
382                k2 = child % childarray(3,2,2,nb,ndir)
383                l1 = child % childarray(4,1,2,nb,ndir)
384                l2 = child % childarray(4,2,2,nb,ndir)
385
386                call procname(parray_temp(kindex_2d(ndir,nb)),i1,i2,j1,j2,k1,k2,l1,l2,.FALSE.,coords(nb),ndir)
387
388            case(5)
389                i1 = child % childarray(1,1,2,nb,ndir)
390                i2 = child % childarray(1,2,2,nb,ndir)
391                j1 = child % childarray(2,1,2,nb,ndir)
392                j2 = child % childarray(2,2,2,nb,ndir)
393                k1 = child % childarray(3,1,2,nb,ndir)
394                k2 = child % childarray(3,2,2,nb,ndir)
395                l1 = child % childarray(4,1,2,nb,ndir)
396                l2 = child % childarray(4,2,2,nb,ndir)
397                m1 = child % childarray(5,1,2,nb,ndir)
398                m2 = child % childarray(5,2,2,nb,ndir)
399
400                call procname(parray5(i1:i2,j1:j2,k1:k2,l1:l2,m1:m2),       &
401                              i1,i2,j1,j2,k1,k2,l1,l2,m1,m2, .FALSE.,coords(nb),ndir)
402            case(6)
403                i1 = child % childarray(1,1,2,nb,ndir)
404                i2 = child % childarray(1,2,2,nb,ndir)
405                j1 = child % childarray(2,1,2,nb,ndir)
406                j2 = child % childarray(2,2,2,nb,ndir)
407                k1 = child % childarray(3,1,2,nb,ndir)
408                k2 = child % childarray(3,2,2,nb,ndir)
409                l1 = child % childarray(4,1,2,nb,ndir)
410                l2 = child % childarray(4,2,2,nb,ndir)
411                m1 = child % childarray(5,1,2,nb,ndir)
412                m2 = child % childarray(5,2,2,nb,ndir)
413                n1 = child % childarray(6,1,2,nb,ndir)
414                n2 = child % childarray(6,2,2,nb,ndir)
415
416                call procname(parray6(i1:i2,j1:j2,k1:k2,l1:l2,m1:m2,n1:n2), &
417                              i1,i2,j1,j2,k1,k2,l1,l2,m1,m2,n1,n2, .FALSE.,coords(nb),ndir)
418            end select
419        endif
420    enddo
421    enddo
422
423    else
424
425    do nb = 1,nbdim
426    do ndir = 1,2
427        if ( (loctab_child(nb) /= (-ndir)) .AND. (loctab_child(nb) /= -3) .AND. child%memberin(nb,ndir) ) then
428            select case(nbdim)
429            case(1)
430                i1 = child % childarray(1,1,2,nb,ndir)
431                i2 = child % childarray(1,2,2,nb,ndir)
432
433                call procname(parray1(i1:i2),                               &
434                              i1,i2, .FALSE.,coords(nb),ndir)
435            case(2)
436                i1 = child % childarray(1,1,2,nb,ndir)
437                i2 = child % childarray(1,2,2,nb,ndir)
438                j1 = child % childarray(2,1,2,nb,ndir)
439                j2 = child % childarray(2,2,2,nb,ndir)
440
441                call procname(parray2(i1:i2,j1:j2),                         &
442                              i1,i2,j1,j2, .FALSE.,coords(nb),ndir)
443            case(3)
444
445                i1 = child % childarray(1,1,2,nb,ndir)
446                i2 = child % childarray(1,2,2,nb,ndir)
447                j1 = child % childarray(2,1,2,nb,ndir)
448                j2 = child % childarray(2,2,2,nb,ndir)
449                k1 = child % childarray(3,1,2,nb,ndir)
450                k2 = child % childarray(3,2,2,nb,ndir)
451
452                call procname(parray3(i1:i2,j1:j2,k1:k2),                   &
453                              i1,i2,j1,j2,k1,k2, .FALSE.,coords(nb),ndir)
454
455            case(4)
456                i1 = child % childarray(1,1,2,nb,ndir)
457                i2 = child % childarray(1,2,2,nb,ndir)
458                j1 = child % childarray(2,1,2,nb,ndir)
459                j2 = child % childarray(2,2,2,nb,ndir)
460                k1 = child % childarray(3,1,2,nb,ndir)
461                k2 = child % childarray(3,2,2,nb,ndir)
462                l1 = child % childarray(4,1,2,nb,ndir)
463                l2 = child % childarray(4,2,2,nb,ndir)
464
465                call procname(parray4(i1:i2,j1:j2,k1:k2,l1:l2),             &
466                              i1,i2,j1,j2,k1,k2,l1,l2, .FALSE.,coords(nb),ndir)
467
468            case(5)
469                i1 = child % childarray(1,1,2,nb,ndir)
470                i2 = child % childarray(1,2,2,nb,ndir)
471                j1 = child % childarray(2,1,2,nb,ndir)
472                j2 = child % childarray(2,2,2,nb,ndir)
473                k1 = child % childarray(3,1,2,nb,ndir)
474                k2 = child % childarray(3,2,2,nb,ndir)
475                l1 = child % childarray(4,1,2,nb,ndir)
476                l2 = child % childarray(4,2,2,nb,ndir)
477                m1 = child % childarray(5,1,2,nb,ndir)
478                m2 = child % childarray(5,2,2,nb,ndir)
479
480                call procname(parray5(i1:i2,j1:j2,k1:k2,l1:l2,m1:m2),       &
481                              i1,i2,j1,j2,k1,k2,l1,l2,m1,m2, .FALSE.,coords(nb),ndir)
482            case(6)
483                i1 = child % childarray(1,1,2,nb,ndir)
484                i2 = child % childarray(1,2,2,nb,ndir)
485                j1 = child % childarray(2,1,2,nb,ndir)
486                j2 = child % childarray(2,2,2,nb,ndir)
487                k1 = child % childarray(3,1,2,nb,ndir)
488                k2 = child % childarray(3,2,2,nb,ndir)
489                l1 = child % childarray(4,1,2,nb,ndir)
490                l2 = child % childarray(4,2,2,nb,ndir)
491                m1 = child % childarray(5,1,2,nb,ndir)
492                m2 = child % childarray(5,2,2,nb,ndir)
493                n1 = child % childarray(6,1,2,nb,ndir)
494                n2 = child % childarray(6,2,2,nb,ndir)
495
496                call procname(parray6(i1:i2,j1:j2,k1:k2,l1:l2,m1:m2,n1:n2), &
497                              i1,i2,j1,j2,k1,k2,l1,l2,m1,m2,n1,n2, .FALSE.,coords(nb),ndir)
498            end select
499        endif
500    enddo
501    enddo
502
503    endif
504!
505
506!---------------------------------------------------------------------------------------------------
507end subroutine Agrif_Correctnd
508!===================================================================================================
509!
510!===================================================================================================
511!  subroutine saveAfterInterp
512!
513!> saves the values of the grid variable on the fine grid after the space interpolation
514!---------------------------------------------------------------------------------------------------
515subroutine saveAfterInterp ( child_var, bounds, kindex, newsize, nbdim )
516!---------------------------------------------------------------------------------------------------
517    TYPE (Agrif_Variable),       INTENT(inout)  :: child_var !< The fine grid variable
518    INTEGER, DIMENSION(nbdim,2), INTENT(in)     :: bounds
519    INTEGER,                     INTENT(inout)  :: kindex    !< Index indicating where this safeguard
520                                                             !<     is done on the fine grid
521    INTEGER,                     INTENT(in)     :: newsize
522    INTEGER,                     INTENT(in)     :: nbdim
523!
524    INTEGER     :: ir,jr,kr,lr,mr,nr
525!
526!    Allocation of the array oldvalues2d
527!
528    if (newsize .LE. 0) return
529!
530    Call Agrif_Checksize(child_var,kindex+newsize)
531
532    if (child_var % interpIndex /= Agrif_Curgrid % parent % ngridstep ) then
533        child_var % oldvalues2d(1,kindex:kindex+newsize-1) = &
534        child_var % oldvalues2d(2,kindex:kindex+newsize-1)
535    endif
536
537    SELECT CASE (nbdim)
538    CASE (1)
539!CDIR ALTCODE
540        do ir = bounds(1,1), bounds(1,2)
541            child_var % oldvalues2d(2,kindex) = parray1(ir)
542            kindex = kindex + 1
543        enddo
544!
545    CASE (2)
546        do jr = bounds(2,1),bounds(2,2)
547!CDIR ALTCODE
548        do ir = bounds(1,1),bounds(1,2)
549            child_var % oldvalues2d(2,kindex) = parray2(ir,jr)
550            kindex = kindex + 1
551        enddo
552        enddo
553!
554    CASE (3)
555        do kr = bounds(3,1),bounds(3,2)
556        do jr = bounds(2,1),bounds(2,2)
557!CDIR ALTCODE
558        do ir = bounds(1,1),bounds(1,2)
559            child_var % oldvalues2d(2,kindex) = parray3(ir,jr,kr)
560            kindex = kindex + 1
561        enddo
562        enddo
563        enddo
564!
565      CASE (4)
566        do lr = bounds(4,1),bounds(4,2)
567        do kr = bounds(3,1),bounds(3,2)
568        do jr = bounds(2,1),bounds(2,2)
569!CDIR ALTCODE
570        do ir = bounds(1,1),bounds(1,2)
571            child_var % oldvalues2d(2,kindex) = parray4(ir,jr,kr,lr)
572            kindex = kindex + 1
573        enddo
574        enddo
575        enddo
576        enddo
577!
578    CASE (5)
579        do mr = bounds(5,1),bounds(5,2)
580        do lr = bounds(4,1),bounds(4,2)
581        do kr = bounds(3,1),bounds(3,2)
582        do jr = bounds(2,1),bounds(2,2)
583!CDIR ALTCODE
584        do ir = bounds(1,1),bounds(1,2)
585            child_var % oldvalues2d(2,kindex) = parray5(ir,jr,kr,lr,mr)
586            kindex = kindex + 1
587        enddo
588        enddo
589        enddo
590        enddo
591        enddo
592!
593    CASE (6)
594        do nr = bounds(6,1),bounds(6,2)
595        do mr = bounds(5,1),bounds(5,2)
596        do lr = bounds(4,1),bounds(4,2)
597        do kr = bounds(3,1),bounds(3,2)
598        do jr = bounds(2,1),bounds(2,2)
599!CDIR ALTCODE
600        do ir = bounds(1,1),bounds(1,2)
601            child_var % oldvalues2d(2,kindex) = parray6(ir,jr,kr,lr,mr,nr)
602            kindex = kindex + 1
603        enddo
604        enddo
605        enddo
606        enddo
607        enddo
608        enddo
609    END SELECT
610!---------------------------------------------------------------------------------------------------
611end subroutine saveAfterInterp
612!===================================================================================================
613!
614!===================================================================================================
615!  subroutine timeInterpolation
616!
617!> subroutine for a linear time interpolation on the child grid
618!---------------------------------------------------------------------------------------------------
619subroutine timeInterpolation ( child_var, bounds, kindex, c1t, c2t, nbdim )
620!---------------------------------------------------------------------------------------------------
621    TYPE (Agrif_Variable)       :: child_var !< The fine grid variable
622    INTEGER, DIMENSION(nbdim,2) :: bounds
623    INTEGER                     :: kindex    !< Index indicating the values of the fine grid got
624                                             !<    before and after the space interpolation and
625                                             !<    used for the time interpolation
626    REAL                        :: c1t, c2t  !< Coefficients for the time interpolation (c2t=1-c1t)
627    INTEGER                     :: nbdim
628!
629    INTEGER :: ir,jr,kr,lr,mr,nr
630    INTEGER :: kindexmax, isize,i
631    REAL,DIMENSION(:),ALLOCATABLE :: tabtemp
632
633    isize = 1
634    DO i=1,nbdim
635      isize = isize * (bounds(i,2)-bounds(i,1)+1)
636    ENDDO
637    IF (isize <= 0) RETURN
638
639    kindexmax = kindex + isize - 1
640    IF (.NOT.ALLOCATED(parray_temp)) THEN
641      ALLOCATE(parray_temp(kindexmax)) 
642    ELSE
643      IF (size(parray_temp) < kindexmax) THEN
644         ALLOCATE(tabtemp(size(parray_temp)))
645         tabtemp = parray_temp
646         DEALLOCATE(parray_temp)
647         ALLOCATE(parray_temp(kindexmax))
648         parray_temp(1:size(tabtemp)) = tabtemp
649         DEALLOCATE(tabtemp)
650      ENDIF
651    ENDIF
652 
653!
654    SELECT CASE (nbdim)
655    CASE (1)
656!CDIR ALTCODE
657        do ir = bounds(1,1),bounds(1,2)
658            parray1(ir) = c2t*child_var % oldvalues2d(1,kindex) + &
659                                      c1t*child_var % oldvalues2d(2,kindex)
660            kindex = kindex + 1
661        enddo
662!
663    CASE (2)
664        do jr = bounds(2,1),bounds(2,2)
665!CDIR ALTCODE
666        do ir = bounds(1,1),bounds(1,2)
667            parray2(ir,jr) = c2t*child_var % oldvalues2d(1,kindex) + &
668                                         c1t*child_var % oldvalues2d(2,kindex)
669            kindex = kindex + 1
670        enddo
671        enddo
672!
673    CASE (3)
674
675        parray_temp(kindex:kindexmax) = c2t*child_var % oldvalues2d(1,kindex:kindexmax) + &
676                                        c1t*child_var % oldvalues2d(2,kindex:kindexmax)
677
678!
679    CASE (4)
680
681        parray_temp(kindex:kindexmax) = c2t*child_var % oldvalues2d(1,kindex:kindexmax) + &
682                                        c1t*child_var % oldvalues2d(2,kindex:kindexmax)
683
684!
685    CASE (5)
686        do mr=bounds(5,1),bounds(5,2)
687        do lr=bounds(4,1),bounds(4,2)
688        do kr=bounds(3,1),bounds(3,2)
689        do jr=bounds(2,1),bounds(2,2)
690!CDIR ALTCODE
691        do ir=bounds(1,1),bounds(1,2)
692            parray5(ir,jr,kr,lr,mr) = c2t*child_var % oldvalues2d(1,kindex) + &
693                                                  c1t*child_var % oldvalues2d(2,kindex)
694            kindex = kindex + 1
695        enddo
696        enddo
697        enddo
698        enddo
699        enddo
700!
701    CASE (6)
702        do nr=bounds(6,1),bounds(6,2)
703        do mr=bounds(5,1),bounds(5,2)
704        do lr=bounds(4,1),bounds(4,2)
705        do kr=bounds(3,1),bounds(3,2)
706        do jr=bounds(2,1),bounds(2,2)
707!CDIR ALTCODE
708        do ir=bounds(1,1),bounds(1,2)
709            parray6(ir,jr,kr,lr,mr,nr) = c2t*child_var % oldvalues2d(1,kindex) + &
710                                                     c1t*child_var % oldvalues2d(2,kindex)
711            kindex = kindex + 1
712        enddo
713        enddo
714        enddo
715        enddo
716        enddo
717        enddo
718    END SELECT
719
720    kindex = kindexmax + 1
721
722!---------------------------------------------------------------------------------------------------
723end subroutine timeInterpolation
724!===================================================================================================
725!
726!===================================================================================================
727!  subroutine Agrif_Checksize
728!
729!> subroutine used in the saveAfterInterp procedure to allocate the oldvalues2d array
730!---------------------------------------------------------------------------------------------------
731subroutine Agrif_Checksize ( child_var, newsize )
732!---------------------------------------------------------------------------------------------------
733    TYPE (Agrif_Variable), INTENT(inout)  :: child_var !< The fine grid variable
734    INTEGER              , INTENT(in)     :: newsize   !< Size of the domains where the boundary
735                                                       !<    conditions are calculated
736!
737    REAL, DIMENSION(:,:), Allocatable :: tempoldvalues ! Temporary array
738!
739    if (.NOT. associated(child_var % oldvalues2d)) then
740!
741        allocate(child_var % oldvalues2d(2,newsize))
742        child_var % oldvalues2d = 0.
743!
744    else
745!
746        if (SIZE(child_var % oldvalues2d,2) < newsize) then
747!
748            allocate(tempoldvalues(2,SIZE(child_var % oldvalues2d,2)))
749            tempoldvalues = child_var % oldvalues2d
750            deallocate(child_var % oldvalues2d)
751            allocate(  child_var % oldvalues2d(2,newsize))
752            child_var % oldvalues2d = 0.
753            child_var % oldvalues2d(:,1:SIZE(tempoldvalues,2)) = tempoldvalues(:,:)
754            deallocate(tempoldvalues)
755!
756        endif
757!
758    endif
759!---------------------------------------------------------------------------------------------------
760end subroutine Agrif_Checksize
761!===================================================================================================
762!
763end module Agrif_Boundary
764
Note: See TracBrowser for help on using the repository browser.