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

source: vendors/AGRIF/CMEMS_2020/AGRIF_FILES/modbcfunction.F90 @ 10725

Last change on this file since 10725 was 10725, checked in by rblod, 5 years ago

Update agrif library and conv see ticket #2129

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