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/dev/AGRIF_FILES – NEMO

source: vendors/AGRIF/dev/AGRIF_FILES/modbcfunction.F90

Last change on this file was 14975, checked in by jchanut, 3 years ago

#2638, merge new AGRIF library into trunk

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