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 vendors/AGRIF/release-4.0.1/AGRIF_FILES – NEMO

source: vendors/AGRIF/release-4.0.1/AGRIF_FILES/modbcfunction.F90 @ 13680

Last change on this file since 13680 was 5656, checked in by timgraham, 9 years ago

Merge of AGRIF branch (branches/2014/dev_r4765_CNRS_agrif) onto the trunk

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