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.
modbcfunction.F90 in branches/2017/dev_r8126_UKMO_AGRIF_vert_interp/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/2017/dev_r8126_UKMO_AGRIF_vert_interp/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modbcfunction.F90 @ 8138

Last change on this file since 8138 was 8138, checked in by timgraham, 7 years ago

Modifications to AGRIF_FILES as received from Laurent (need to check that these are definitely needed)

  • Property svn:keywords set to Id
File size: 31.5 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!> Module Agrif_BcFunction.
24!
25module Agrif_BcFunction
26!
27!     Modules used:
28!
29    use Agrif_Boundary
30    use Agrif_Update
31    use Agrif_Save
32!
33    implicit none
34!
35    interface Agrif_Set_Parent
36        module procedure Agrif_Set_Parent_int,      &
37                         Agrif_Set_Parent_real4,    &
38                         Agrif_Set_Parent_real8
39    end interface
40!
41    interface Agrif_Save_Forrestore
42        module procedure Agrif_Save_Forrestore0d,   &
43                         Agrif_Save_Forrestore2d,   &
44                         Agrif_Save_Forrestore3d,   &
45                         Agrif_Save_Forrestore4d
46    end interface
47!
48contains
49!
50!===================================================================================================
51!  subroutine Agrif_Set_parent_int
52!
53!> To set the TYPE of the variable
54!---------------------------------------------------------------------------------------------------
55subroutine Agrif_Set_parent_int(tabvarsindic,value)
56!---------------------------------------------------------------------------------------------------
57    integer, intent(in)     :: tabvarsindic !< indice of the variable in tabvars
58    integer, intent(in)     :: value        !< input value
59!
60
61integer :: i
62logical :: i_found
63
64i_found = .FALSE.
65
66do i=1,Agrif_NbVariables(4)
67  if (LOC(tabvarsindic) == LOC(agrif_curgrid%tabvars_i(i)%iarray0)) then
68     agrif_curgrid%tabvars_i(i)%parent_var%iarray0 = value
69     i_found = .TRUE.
70     EXIT
71  endif
72enddo
73
74if (.NOT.i_found) STOP 'Agrif_Set_Integer : Variable not found'
75
76!---------------------------------------------------------------------------------------------------
77end subroutine Agrif_Set_parent_int
78!===================================================================================================
79!
80!===================================================================================================
81!  subroutine Agrif_Set_parent_real4
82!---------------------------------------------------------------------------------------------------
83!> To set the TYPE of the variable
84!---------------------------------------------------------------------------------------------------
85subroutine Agrif_Set_parent_real4 ( tabvarsindic, value )
86!---------------------------------------------------------------------------------------------------
87    real(kind=4), intent(in)     :: tabvarsindic !< input variable
88    real(kind=4),intent(in) :: value        !< input value for the parent grid
89
90integer :: i
91logical :: i_found
92
93i_found = .FALSE.
94
95do i=1,Agrif_NbVariables(2)
96  if (LOC(tabvarsindic) == LOC(agrif_curgrid%tabvars_r(i)%array0)) then
97     agrif_curgrid%tabvars_r(i)%parent_var%array0 = value
98     agrif_curgrid%tabvars_r(i)%parent_var%sarray0 = value
99     i_found = .TRUE.
100     EXIT
101  endif
102enddo
103
104IF (.NOT.i_found) THEN
105do i=1,Agrif_NbVariables(2)
106  if (LOC(tabvarsindic) == LOC(agrif_curgrid%tabvars_r(i)%sarray0)) then
107     agrif_curgrid%tabvars_r(i)%parent_var%array0 = value
108     agrif_curgrid%tabvars_r(i)%parent_var%sarray0 = value
109     i_found = .TRUE.
110     EXIT
111  endif
112enddo
113ENDIF
114
115if (.NOT.i_found) STOP 'Agrif_Set_parent_real4 : Variable not found'
116!---------------------------------------------------------------------------------------------------
117end subroutine Agrif_Set_parent_real4
118!===================================================================================================
119!
120!===================================================================================================
121!  subroutine Agrif_Set_parent_real8
122!---------------------------------------------------------------------------------------------------
123!> To set the TYPE of the variable
124!---------------------------------------------------------------------------------------------------
125subroutine Agrif_Set_parent_real8 ( tabvarsindic, value )
126!---------------------------------------------------------------------------------------------------
127    real(kind=8), intent(in)     :: tabvarsindic !< input variable
128    real(kind=8),intent(in) :: value        !< input value for the parent grid
129
130integer :: i
131logical :: i_found
132
133i_found = .FALSE.
134
135do i=1,Agrif_NbVariables(2)
136  if (LOC(tabvarsindic) == LOC(agrif_curgrid%tabvars_r(i)%array0)) then
137     agrif_curgrid%tabvars_r(i)%parent_var%darray0 = value
138     agrif_curgrid%tabvars_r(i)%parent_var%array0 = value
139     i_found = .TRUE.
140     EXIT
141  endif
142enddo
143
144IF (.NOT.i_found) THEN
145do i=1,Agrif_NbVariables(2)
146  if (LOC(tabvarsindic) == LOC(agrif_curgrid%tabvars_r(i)%darray0)) then
147     agrif_curgrid%tabvars_r(i)%parent_var%darray0 = value
148     agrif_curgrid%tabvars_r(i)%parent_var%array0 = value
149     i_found = .TRUE.
150     EXIT
151  endif
152enddo
153ENDIF
154
155if (.NOT.i_found) STOP 'Agrif_Set_parent_real8 : Variable not found'
156
157!---------------------------------------------------------------------------------------------------
158end subroutine Agrif_Set_parent_real8
159!===================================================================================================
160!
161!===================================================================================================
162!  subroutine Agrif_Set_bc
163!---------------------------------------------------------------------------------------------------
164subroutine Agrif_Set_bc ( tabvarsindic, bcinfsup, Interpolationshouldbemade )
165!---------------------------------------------------------------------------------------------------
166    integer,               intent(in)   :: tabvarsindic !< indice of the variable in tabvars
167    integer, dimension(2), intent(in)   :: bcinfsup     !< bcinfsup
168    logical, optional,     intent(in)   :: Interpolationshouldbemade !< interpolation should be made
169!
170    integer                         :: indic ! indice of the variable in tabvars
171    type(Agrif_Variable),  pointer  :: var
172!
173    indic = Agrif_Curgrid % tabvars_i(tabvarsindic) % iarray0
174!
175    if (indic <= 0) then
176        var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
177    else
178        print*,"Agrif_Set_bc : warning indic >= 0 !!!"
179        var => Agrif_Curgrid % tabvars(indic)
180    endif
181
182    if (.not.associated(var)) return ! Grand mother grid case
183!
184    if ( Agrif_Curgrid % fixedrank /= 0 ) then
185        if ( .not.associated(var % oldvalues2D) ) then
186            allocate(var % oldvalues2D(2,1))
187            var % interpIndex = -1
188            var % oldvalues2D = 0.
189        endif
190        if ( present(Interpolationshouldbemade) ) then
191            var % Interpolationshouldbemade = Interpolationshouldbemade
192        endif
193    endif
194!
195    var % bcinf = bcinfsup(1)
196    var % bcsup = bcinfsup(2)
197!---------------------------------------------------------------------------------------------------
198end subroutine Agrif_Set_bc
199!===================================================================================================
200!
201!===================================================================================================
202!  subroutine Agrif_Set_interp
203!---------------------------------------------------------------------------------------------------
204subroutine Agrif_Set_interp ( tabvarsindic, interp, interp1, interp2, interp3 , interp4)
205!---------------------------------------------------------------------------------------------------
206    integer,           intent(in)   :: tabvarsindic !< indice of the variable in tabvars
207    integer, optional, intent(in)   :: interp, interp1, interp2, interp3, interp4
208!
209    integer                         :: indic ! indice of the variable in tabvars
210    type(Agrif_Variable), pointer   :: var
211!
212    indic = Agrif_Curgrid % tabvars_i(tabvarsindic) % iarray0
213!
214    if (indic <= 0) then
215        var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
216    else
217        print*,"Agrif_Set_interp : warning indic >= 0 !!!"
218        var => Agrif_Mygrid % tabvars(indic)
219    endif
220!
221    var % type_interp = Agrif_Constant
222!
223    if (present(interp))    var % type_interp    = interp
224    if (present(interp1))   var % type_interp(1) = interp1
225    if (present(interp2))   var % type_interp(2) = interp2
226    if (present(interp3))   var % type_interp(3) = interp3
227    if (present(interp4))   var % type_interp(4) = interp4
228!---------------------------------------------------------------------------------------------------
229end subroutine Agrif_Set_interp
230!===================================================================================================
231!
232!===================================================================================================
233!  subroutine Agrif_Set_bcinterp
234!---------------------------------------------------------------------------------------------------
235subroutine Agrif_Set_bcinterp ( tabvarsindic, interp,   interp1,  interp2,  interp3, interp4, &
236                                              interp11, interp12, interp21, interp22 )
237!---------------------------------------------------------------------------------------------------
238    INTEGER,           intent(in)   :: tabvarsindic !< indice of the variable in tabvars
239    INTEGER, OPTIONAL, intent(in)   :: interp,   interp1,  interp2,  interp3, interp4
240    INTEGER, OPTIONAL, intent(in)   :: interp11, interp12, interp21, interp22
241!
242    INTEGER                         :: indic ! indice of the variable in tabvars
243    TYPE(Agrif_Variable), pointer   :: var
244!
245    indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
246!
247    if (indic <= 0) then
248        var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
249    else
250        print*,"Agrif_Set_bcinterp : warning indic >= 0 !!!"
251        var => Agrif_Mygrid % tabvars(indic)
252    endif
253!
254    var % type_interp_bc = Agrif_Constant
255!
256    if (present(interp))    var % type_interp_bc      = interp
257    if (present(interp1))   var % type_interp_bc(:,1) = interp1
258    if (present(interp11))  var % type_interp_bc(1,1) = interp11
259    if (present(interp12))  var % type_interp_bc(1,2) = interp12
260    if (present(interp2))   var % type_interp_bc(:,2) = interp2
261    if (present(interp21))  var % type_interp_bc(2,1) = interp21
262    if (present(interp22))  var % type_interp_bc(2,2) = interp22
263    if (present(interp3))   var % type_interp_bc(:,3) = interp3
264    if (present(interp4))   var % type_interp_bc(:,4) = interp4
265!---------------------------------------------------------------------------------------------------
266end subroutine Agrif_Set_bcinterp
267!===================================================================================================
268!
269!===================================================================================================
270!  subroutine Agrif_Set_UpdateType
271!---------------------------------------------------------------------------------------------------
272subroutine Agrif_Set_UpdateType ( tabvarsindic, update,  update1, update2, &
273                                                update3, update4, update5 )
274!---------------------------------------------------------------------------------------------------
275    INTEGER,           intent(in) :: tabvarsindic     !< indice of the variable in tabvars
276    INTEGER, OPTIONAL, intent(in) :: update, update1, update2, update3, update4, update5
277!
278    INTEGER                         :: indic ! indice of the variable in tabvars
279    type(Agrif_Variable),  pointer  :: root_var
280!
281    indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
282!
283    if (indic <= 0) then
284        root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
285    else
286        print*,"Agrif_Set_UpdateType : warning indic >= 0 !!!"
287        root_var => Agrif_Mygrid % tabvars(indic)
288    endif
289!
290    root_var % type_update = Agrif_Update_Copy
291    if (present(update))    root_var % type_update    = update
292    if (present(update1))   root_var % type_update(1) = update1
293    if (present(update2))   root_var % type_update(2) = update2
294    if (present(update3))   root_var % type_update(3) = update3
295    if (present(update4))   root_var % type_update(4) = update4
296    if (present(update5))   root_var % type_update(5) = update5
297!---------------------------------------------------------------------------------------------------
298end subroutine Agrif_Set_UpdateType
299!===================================================================================================
300!
301!===================================================================================================
302!  subroutine Agrif_Set_restore
303!---------------------------------------------------------------------------------------------------
304subroutine Agrif_Set_restore ( tabvarsindic )
305!---------------------------------------------------------------------------------------------------
306    INTEGER, intent(in) :: tabvarsindic     !< indice of the variable in tabvars
307!
308    INTEGER :: indic  !  indice of the variable in tabvars
309!
310    indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
311!
312    Agrif_Mygrid%tabvars(indic) % restore = .TRUE.
313!---------------------------------------------------------------------------------------------------
314end subroutine Agrif_Set_restore
315!===================================================================================================
316!
317!===================================================================================================
318!  subroutine Agrif_Init_variable
319!---------------------------------------------------------------------------------------------------
320subroutine Agrif_Init_variable ( tabvarsindic, procname )
321!---------------------------------------------------------------------------------------------------
322    INTEGER, intent(in)  :: tabvarsindic     !< indice of the variable in tabvars
323    procedure()          :: procname         !< Data recovery procedure
324!
325    if ( Agrif_Curgrid%level <= 0 ) return
326!
327    call Agrif_Interp_variable(tabvarsindic, procname)
328    call Agrif_Bc_variable(tabvarsindic, procname, 1.)
329!---------------------------------------------------------------------------------------------------
330end subroutine Agrif_Init_variable
331!===================================================================================================
332!
333!===================================================================================================
334!  subroutine Agrif_Bc_variable
335!---------------------------------------------------------------------------------------------------
336subroutine Agrif_Bc_variable ( tabvarsindic, procname, calledweight )
337!---------------------------------------------------------------------------------------------------
338    integer,        intent(in) :: tabvarsindic     !< indice of the variable in tabvars
339    procedure()                :: procname
340    real, optional, intent(in) :: calledweight
341!
342    real    :: weight
343    logical :: pweight
344    integer :: indic
345    integer :: nbdim
346    type(Agrif_Variable), pointer :: root_var
347    type(Agrif_Variable), pointer :: parent_var
348    type(Agrif_Variable), pointer :: child_var
349    type(Agrif_Variable), pointer :: child_tmp      ! Temporary variable on the child grid
350!
351    if ( Agrif_Curgrid%level <= 0 ) return
352!
353    indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
354!
355    if ( present(calledweight) ) then
356        weight  = calledweight
357        pweight = .true.
358    else
359        weight  = 0.
360        pweight = .false.
361    endif
362!
363    if (indic <= 0) then
364        child_var  => Agrif_Search_Variable(Agrif_Curgrid,-indic)
365        parent_var => child_var % parent_var
366        root_var   => child_var % root_var
367    else
368        print*,"Agrif_Bc_variable : warning indic >= 0 !!!"
369        child_var  => Agrif_Curgrid % tabvars(indic)
370        parent_var => Agrif_Curgrid % parent % tabvars(indic)
371        root_var   => Agrif_Mygrid % tabvars(indic)
372    endif
373!
374    nbdim = root_var % nbdim
375!
376    select case( nbdim )
377    case(1)
378        allocate(parray1(child_var%lb(1):child_var%ub(1)))
379    case(2)
380        allocate(parray2(child_var%lb(1):child_var%ub(1), &
381                         child_var%lb(2):child_var%ub(2) ))
382    case(3)
383        allocate(parray3(child_var%lb(1):child_var%ub(1), &
384                         child_var%lb(2):child_var%ub(2), &
385                         child_var%lb(3):child_var%ub(3) ))
386    case(4)
387        allocate(parray4(child_var%lb(1):child_var%ub(1), &
388                         child_var%lb(2):child_var%ub(2), &
389                         child_var%lb(3):child_var%ub(3), &
390                         child_var%lb(4):child_var%ub(4) ))
391    case(5)
392        allocate(parray5(child_var%lb(1):child_var%ub(1), &
393                         child_var%lb(2):child_var%ub(2), &
394                         child_var%lb(3):child_var%ub(3), &
395                         child_var%lb(4):child_var%ub(4), &
396                         child_var%lb(5):child_var%ub(5) ))
397    case(6)
398        allocate(parray6(child_var%lb(1):child_var%ub(1), &
399                         child_var%lb(2):child_var%ub(2), &
400                         child_var%lb(3):child_var%ub(3), &
401                         child_var%lb(4):child_var%ub(4), &
402                         child_var%lb(5):child_var%ub(5), &
403                         child_var%lb(6):child_var%ub(6) ))
404    end select
405!
406!   Create temporary child variable
407    allocate(child_tmp)
408!
409    child_tmp % root_var => root_var
410    child_tmp % oldvalues2D => child_var % oldvalues2D
411!
412!   Index indicating if a space interpolation is necessary
413    child_tmp % interpIndex =  child_var % interpIndex
414    child_tmp % list_interp => child_var % list_interp
415    child_tmp % Interpolationshouldbemade = child_var % Interpolationshouldbemade
416!
417    child_tmp % point = child_var % point
418    child_tmp % lb = child_var % lb
419    child_tmp % ub = child_var % ub
420!
421    child_tmp % bcinf = child_var % bcinf
422    child_tmp % bcsup = child_var % bcsup
423!
424    child_tmp % childarray = child_var % childarray
425    child_tmp % memberin   = child_var % memberin
426!
427    call Agrif_CorrectVariable(parent_var, child_tmp, pweight, weight, procname)
428!
429    child_var % childarray = child_tmp % childarray
430    child_var % memberin   = child_tmp % memberin
431!
432    child_var % oldvalues2D => child_tmp % oldvalues2D
433    child_var % list_interp => child_tmp % list_interp
434!
435    child_var % interpIndex = child_tmp % interpIndex
436!
437    deallocate(child_tmp)
438!
439    select case( nbdim )
440        case(1); deallocate(parray1)
441        case(2); deallocate(parray2)
442        case(3); deallocate(parray3)
443        case(4); deallocate(parray4)
444        case(5); deallocate(parray5)
445        case(6); deallocate(parray6)
446    end select
447!---------------------------------------------------------------------------------------------------
448end subroutine Agrif_Bc_variable
449!===================================================================================================
450!
451!===================================================================================================
452!  subroutine Agrif_Interp_variable
453!---------------------------------------------------------------------------------------------------
454subroutine Agrif_Interp_variable ( tabvarsindic, procname )
455!---------------------------------------------------------------------------------------------------
456    integer,     intent(in)     :: tabvarsindic     !< indice of the variable in tabvars
457    procedure()                 :: procname         !< Data recovery procedure
458!
459    integer :: nbdim
460    integer :: indic  ! indice of the variable in tabvars
461    logical :: torestore
462    type(Agrif_Variable), pointer   :: root_var
463    type(Agrif_Variable), pointer   :: parent_var       ! Variable on the parent grid
464    type(Agrif_Variable), pointer   :: child_var        ! Variable on the parent grid
465    type(Agrif_Variable), pointer   :: child_tmp        ! Temporary variable on the child grid
466!
467    if ( Agrif_Curgrid%level <= 0 ) return
468!
469    indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
470!
471    if (indic <= 0) then
472        child_var  => Agrif_Search_Variable(Agrif_Curgrid,-indic)
473        parent_var => child_var % parent_var
474        root_var   => child_var % root_var
475    else
476        print*,"Agrif_Interp_variable : warning indic >= 0 !!!"
477        child_var  => Agrif_Curgrid % tabvars(indic)
478        parent_var => Agrif_Curgrid % parent % tabvars(indic)
479        root_var   => Agrif_Mygrid % tabvars(indic)
480    endif
481!
482    nbdim     = root_var % nbdim
483    torestore = root_var % restore
484!
485    allocate(child_tmp)
486!
487    child_tmp % root_var => root_var
488    child_tmp % nbdim = root_var % nbdim
489    child_tmp % point = child_var % point
490    child_tmp % lb = child_var % lb
491    child_tmp % ub = child_var % ub
492    child_tmp % interpIndex =  child_var % interpIndex
493    child_tmp % list_interp => child_var % list_interp
494    child_tmp % Interpolationshouldbemade = child_var % Interpolationshouldbemade
495!
496    if ( torestore ) then
497        select case( nbdim )
498        case(1)
499            parray1 = child_var % array1
500            child_tmp % restore1D => child_var % restore1D
501        case(2)
502            parray2 = child_var % array2
503            child_tmp % restore2D => child_var % restore2D
504        case(3)
505            parray3 = child_var % array3
506            child_tmp % restore3D => child_var % restore3D
507        case(4)
508            parray4 = child_var % array4
509            child_tmp % restore4D => child_var % restore4D
510        case(5)
511            parray5 = child_var % array5
512            child_tmp % restore5D => child_var % restore5D
513        case(6)
514            parray6 = child_var % array6
515            child_tmp % restore6D => child_var % restore6D
516        end select
517    endif
518!
519    call Agrif_InterpVariable(parent_var, child_tmp, torestore, procname)
520!
521    child_var % list_interp => child_tmp % list_interp
522!
523    deallocate(child_tmp)
524!---------------------------------------------------------------------------------------------------
525end subroutine Agrif_Interp_variable
526!===================================================================================================
527!
528!===================================================================================================
529!  subroutine Agrif_Update_Variable
530!---------------------------------------------------------------------------------------------------
531subroutine Agrif_Update_Variable ( tabvarsindic, procname, &
532                                   locupdate, locupdate1, locupdate2, locupdate3, locupdate4 )
533!---------------------------------------------------------------------------------------------------
534    integer,               intent(in)           :: tabvarsindic     !< Indice of the variable in tabvars
535    procedure()                                 :: procname         !< Data recovery procedure
536    integer, dimension(2), intent(in), optional :: locupdate
537    integer, dimension(2), intent(in), optional :: locupdate1
538    integer, dimension(2), intent(in), optional :: locupdate2
539    integer, dimension(2), intent(in), optional :: locupdate3
540    integer, dimension(2), intent(in), optional :: locupdate4
541!---------------------------------------------------------------------------------------------------
542    integer :: indic
543    integer :: nbdim
544    integer, dimension(6)           :: updateinf    ! First positions where interpolations are calculated
545    integer, dimension(6)           :: updatesup    ! Last  positions where interpolations are calculated
546    type(Agrif_Variable), pointer   :: root_var
547    type(Agrif_Variable), pointer   :: parent_var
548    type(Agrif_Variable), pointer   :: child_var
549!
550    if ( Agrif_Root() .AND. (.not.agrif_coarse) ) return
551    if (agrif_curgrid%grand_mother_grid) return
552!
553    indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
554!
555    if (indic <= 0) then
556        child_var  => Agrif_Search_Variable(Agrif_Curgrid, -indic)
557        parent_var => child_var % parent_var
558
559        if (.not.associated(parent_var)) then
560          ! can occur during the first update of Agrif_Coarsegrid (if any)
561          parent_var => Agrif_Search_Variable(Agrif_Curgrid % parent, -indic)
562          child_var % parent_var => parent_var
563        endif
564
565        root_var   => child_var % root_var
566    else
567        print*,"Agrif_Update_Variable : warning indic >= 0 !!!"
568        root_var   => Agrif_Mygrid  % tabvars(indic)
569        child_var  => Agrif_Curgrid % tabvars(indic)
570        parent_var => Agrif_Curgrid % parent % tabvars(indic)
571    endif
572!
573    nbdim = root_var % nbdim
574!
575    updateinf = -99
576    updatesup = -99
577!
578    if ( present(locupdate) ) then
579        updateinf(1:nbdim) = locupdate(1)
580        updatesup(1:nbdim) = locupdate(2)
581    endif
582!
583    if ( present(locupdate1) ) then
584        updateinf(1) = locupdate1(1)
585        updatesup(1) = locupdate1(2)
586    endif
587!
588    if ( present(locupdate2) ) then
589        updateinf(2) = locupdate2(1)
590        updatesup(2) = locupdate2(2)
591    endif
592
593    if ( present(locupdate3) ) then
594        updateinf(3) = locupdate3(1)
595        updatesup(3) = locupdate3(2)
596    endif
597
598    if ( present(locupdate4) ) then
599        updateinf(4) = locupdate4(1)
600        updatesup(4) = locupdate4(2)
601    endif
602!
603    call Agrif_UpdateVariable( parent_var, child_var, updateinf, updatesup, procname )
604!---------------------------------------------------------------------------------------------------
605end subroutine Agrif_Update_Variable
606!===================================================================================================
607!
608!===================================================================================================
609!  subroutine Agrif_Save_ForRestore0D
610!---------------------------------------------------------------------------------------------------
611subroutine Agrif_Save_ForRestore0D ( tabvarsindic0, tabvarsindic )
612!---------------------------------------------------------------------------------------------------
613    integer, intent(in) :: tabvarsindic0, tabvarsindic
614!
615    type(Agrif_Variable), pointer   :: root_var, save_var
616    integer :: nbdim
617!
618    root_var => Agrif_Mygrid % tabvars(tabvarsindic0)
619    save_var => Agrif_Curgrid % tabvars(tabvarsindic0)
620    nbdim =  root_var % nbdim
621!
622    select case(nbdim)
623        case(2); call Agrif_Save_ForRestore2D(save_var % array2, tabvarsindic)
624        case(3); call Agrif_Save_ForRestore3D(save_var % array3, tabvarsindic)
625        case(4); call Agrif_Save_ForRestore4D(save_var % array4, tabvarsindic)
626    end select
627!---------------------------------------------------------------------------------------------------
628end subroutine Agrif_Save_ForRestore0D
629!===================================================================================================
630!
631!===================================================================================================
632!  subroutine Agrif_Save_ForRestore2D
633!---------------------------------------------------------------------------------------------------
634subroutine Agrif_Save_ForRestore2D ( q, tabvarsindic )
635!---------------------------------------------------------------------------------------------------
636    real, dimension(:,:), intent(in)        :: q
637    integer,              intent(in)        :: tabvarsindic
638!
639    type(Agrif_Variable),  pointer  :: root_var, save_var
640    integer                         :: indic
641!
642    indic = tabvarsindic
643    if (tabvarsindic >= 0) then
644        if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then
645            indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
646        endif
647    endif
648!
649    if (indic <= 0) then
650        save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
651        root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
652    else
653        save_var => Agrif_Curgrid % tabvars(indic)
654        root_var => Agrif_Mygrid % tabvars(indic)
655    endif
656!
657    if ( .not.allocated(save_var%array2) ) then
658        allocate(save_var%array2(save_var%lb(1):save_var%ub(1),  &
659                                 save_var%lb(2):save_var%ub(2)))
660    endif
661!
662    save_var % array2  = q
663    root_var % restore = .true.
664!---------------------------------------------------------------------------------------------------
665end subroutine Agrif_Save_ForRestore2D
666!===================================================================================================
667!
668!===================================================================================================
669!  subroutine Agrif_Save_ForRestore3D
670!---------------------------------------------------------------------------------------------------
671subroutine Agrif_Save_ForRestore3D ( q, tabvarsindic )
672!---------------------------------------------------------------------------------------------------
673    real, dimension(:,:,:), intent(in)      :: q
674    integer,                intent(in)      :: tabvarsindic
675!
676    type(Agrif_Variable),  pointer  :: root_var, save_var
677    integer                         :: indic
678!
679    indic = tabvarsindic
680    if (tabvarsindic >= 0) then
681        if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then
682            indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
683        endif
684    endif
685!
686    if (indic <= 0) then
687        save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
688        root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
689    else
690        save_var => Agrif_Curgrid % tabvars(indic)
691        root_var => Agrif_Mygrid % tabvars(indic)
692    endif
693!
694    if ( .not.allocated(save_var%array3) ) then
695        allocate(save_var%array3(save_var%lb(1):save_var%ub(1), &
696                                 save_var%lb(2):save_var%ub(2), &
697                                 save_var%lb(3):save_var%ub(3)))
698    endif
699!
700    save_var % array3  = q
701    root_var % restore = .true.
702!---------------------------------------------------------------------------------------------------
703end subroutine Agrif_Save_ForRestore3D
704!===================================================================================================
705!
706!===================================================================================================
707!  subroutine Agrif_Save_ForRestore4D
708!---------------------------------------------------------------------------------------------------
709subroutine Agrif_Save_ForRestore4D ( q, tabvarsindic )
710!---------------------------------------------------------------------------------------------------
711    real, dimension(:,:,:,:), intent(in)    :: q
712    integer,                  intent(in)    :: tabvarsindic
713!
714    type(Agrif_Variable),  pointer  :: root_var, save_var
715    integer                         :: indic
716!
717    indic = tabvarsindic
718    if (tabvarsindic >= 0) then
719        if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then
720            indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
721        endif
722    endif
723!
724    if (indic <= 0) then
725        save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
726        root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
727    else
728        save_var => Agrif_Curgrid % tabvars(indic)
729        root_var => Agrif_Mygrid % tabvars(indic)
730    endif
731!
732    if (.not.allocated(save_var%array4)) then
733        allocate(save_var%array4(save_var%lb(1):save_var%ub(1),&
734                                 save_var%lb(2):save_var%ub(2),&
735                                 save_var%lb(3):save_var%ub(3),&
736                                 save_var%lb(4):save_var%ub(4)))
737    endif
738!
739    save_var % array4  = q
740    root_var % restore = .true.
741!---------------------------------------------------------------------------------------------------
742end subroutine Agrif_Save_ForRestore4D
743!===================================================================================================
744!
745end module Agrif_BcFunction
Note: See TracBrowser for help on using the repository browser.