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

source: vendors/AGRIF/CMEMS_2020/AGRIF_FILES/modutil.F90 @ 10084

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

Compilation for NESTING_AGRIF

  • Property svn:keywords set to Id
File size: 41.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_Util
24!!
25!! This module contains the two procedures called in the main program :
26!! - #Agrif_Init_Grids allows the initialization of the root coarse grid
27!! - #Agrif_Step allows the creation of the grid hierarchy and the management of the time integration.
28!
29module Agrif_Util
30!
31    use Agrif_Clustering
32    use Agrif_BcFunction
33    use Agrif_seq
34    use Agrif_Link
35!
36    implicit none
37!
38    abstract interface
39        subroutine step_proc()
40        end subroutine step_proc
41    end interface
42!
43contains
44!
45!===================================================================================================
46!  subroutine Agrif_Step
47!
48!> Creates the grid hierarchy and manages the time integration procedure.
49!> It is called in the main program.
50!> Calls subroutines #Agrif_Regrid and #Agrif_Integrate.
51!---------------------------------------------------------------------------------------------------
52subroutine Agrif_Step ( procname )
53!---------------------------------------------------------------------------------------------------
54    procedure(step_proc)    :: procname     !< subroutine to call on each grid
55    type(agrif_grid), pointer :: ref_grid
56!
57! Set the clustering variables
58    call Agrif_clustering_def()
59!
60! Creation and initialization of the grid hierarchy
61    if ( Agrif_USE_ONLY_FIXED_GRIDS == 1 ) then
62!
63        if ( (Agrif_Mygrid % ngridstep == 0) .AND. (.not. Agrif_regrid_has_been_done) ) then
64            call Agrif_Regrid()
65            Agrif_regrid_has_been_done = .TRUE.
66        endif
67!
68    else
69!
70        if (mod(Agrif_Mygrid % ngridstep,Agrif_Regridding) == 0) then
71            call Agrif_Regrid()
72        endif
73!
74    endif
75!
76! Time integration of the grid hierarchy
77    if (agrif_coarse) then
78      ref_grid => agrif_coarsegrid
79    else
80      ref_grid => agrif_mygrid
81    endif
82    if ( Agrif_Parallel_sisters ) then
83        call Agrif_Integrate_Parallel(ref_grid,procname)
84    else
85        call Agrif_Integrate(ref_grid,procname)
86    endif
87!
88    if ( ref_grid%child_list%nitems > 0 ) call Agrif_Instance(ref_grid)
89!---------------------------------------------------------------------------------------------------
90end subroutine Agrif_Step
91!===================================================================================================
92!
93!===================================================================================================
94!  subroutine Agrif_Step_Child
95!
96!> Apply 'procname' to each grid of the hierarchy
97!---------------------------------------------------------------------------------------------------
98subroutine Agrif_Step_Child ( procname )
99!---------------------------------------------------------------------------------------------------
100    procedure(step_proc)    :: procname     !< subroutine to call on each grid
101!
102    if ( Agrif_Parallel_sisters ) then
103        call Agrif_Integrate_Child_Parallel(Agrif_Mygrid,procname)
104    else
105        call Agrif_Integrate_Child(Agrif_Mygrid,procname)
106    endif
107!
108    if ( Agrif_Mygrid%child_list%nitems > 0 ) call Agrif_Instance(Agrif_Mygrid)
109!---------------------------------------------------------------------------------------------------
110end subroutine Agrif_Step_Child
111!===================================================================================================
112!
113!===================================================================================================
114!  subroutine Agrif_Regrid
115!
116!> Creates the grid hierarchy from fixed grids and adaptive mesh refinement.
117!---------------------------------------------------------------------------------------------------
118subroutine Agrif_Regrid ( procname )
119!---------------------------------------------------------------------------------------------------
120    procedure(init_proc), optional    :: procname     !< Initialisation subroutine (Default: Agrif_InitValues)
121!
122    type(Agrif_Rectangle), pointer     :: coarsegrid_fixed
123    type(Agrif_Rectangle), pointer     :: coarsegrid_moving
124    integer                            :: i, j
125    integer :: nunit
126    logical :: BEXIST
127    TYPE(Agrif_Rectangle)            :: newrect    ! Pointer on a new grid
128    integer :: is_coarse, rhox, rhoy, rhoz, rhot
129!
130    if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) &
131        call Agrif_detect_all(Agrif_Mygrid)     ! Detection of areas to be refined
132!
133    allocate(coarsegrid_fixed)
134    allocate(coarsegrid_moving)
135!
136    if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) &
137        call Agrif_Cluster_All(Agrif_Mygrid,coarsegrid_moving) ! Clustering
138!
139    if ( Agrif_USE_FIXED_GRIDS == 1 .OR. Agrif_USE_ONLY_FIXED_GRIDS == 1 ) then
140!
141        if (Agrif_Mygrid % ngridstep == 0) then
142!           
143            nunit = Agrif_Get_Unit()
144            open(nunit, file='AGRIF_FixedGrids.in', form='formatted', status="old", ERR=99)
145            if (agrif_coarse) then ! SKIP the coarse grid declaration
146                if (Agrif_Probdim == 3) then
147                    read(nunit,*) is_coarse, rhox, rhoy, rhoz, rhot
148                elseif (Agrif_Probdim == 2) then
149                    read(nunit,*) is_coarse, rhox, rhoy, rhot
150                elseif (Agrif_Probdim == 2) then
151                    read(nunit,*) is_coarse, rhox, rhot
152                endif
153            endif
154!           Creation of the grid hierarchy from the Agrif_FixedGrids.in file
155            do i = 1,Agrif_Probdim
156                coarsegrid_fixed % imin(i) = 1
157                coarsegrid_fixed % imax(i) = Agrif_Mygrid % nb(i) + 1
158            enddo
159            j = 1
160            call Agrif_Read_Fix_Grd(coarsegrid_fixed,j,nunit)
161            close(nunit)
162!
163            call Agrif_gl_clear(Agrif_oldmygrid)
164!
165!           Creation of the grid hierarchy from coarsegrid_fixed
166            call Agrif_Create_Grids(Agrif_Mygrid,coarsegrid_fixed)
167           
168        else
169            call Agrif_gl_copy(Agrif_oldmygrid, Agrif_Mygrid % child_list)
170        endif
171    else
172        call Agrif_gl_copy(Agrif_oldmygrid, Agrif_Mygrid % child_list)
173        call Agrif_gl_clear(Agrif_Mygrid % child_list)
174    endif
175!
176    if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) then
177!
178        call Agrif_Save_All(Agrif_oldmygrid)
179        call Agrif_Free_before_All(Agrif_oldmygrid)
180!
181!       Creation of the grid hierarchy from coarsegrid_moving
182        call Agrif_Create_Grids(Agrif_Mygrid,coarsegrid_moving)
183!
184    endif
185!
186!   Initialization of the grid hierarchy by copy or interpolation
187!
188#if defined AGRIF_MPI
189    if ( Agrif_Parallel_sisters ) then
190        call Agrif_Init_Hierarchy_Parallel_1(Agrif_Mygrid)
191        call Agrif_Init_Hierarchy_Parallel_2(Agrif_Mygrid,procname)
192    else
193        call Agrif_Init_Hierarchy(Agrif_Mygrid,procname)
194    endif
195#else
196        call Agrif_Init_Hierarchy(Agrif_Mygrid,procname)
197#endif
198!
199    if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) call Agrif_Free_after_All(Agrif_oldmygrid)
200!
201    Agrif_regrid_has_been_done = .TRUE.
202!
203    call Agrif_Instance( Agrif_Mygrid )
204!
205    deallocate(coarsegrid_fixed)
206    deallocate(coarsegrid_moving)
207!
208    return
209!
210!     Opening error
211!
21299  INQUIRE(FILE='AGRIF_FixedGrids.in',EXIST=BEXIST)
213    if (.not. BEXIST) then
214        print*,'ERROR : File AGRIF_FixedGrids.in not found.'
215        STOP
216    else
217        print*,'Error opening file AGRIF_FixedGrids.in'
218        STOP
219    endif
220!---------------------------------------------------------------------------------------------------
221end subroutine Agrif_Regrid
222!===================================================================================================
223!
224!===================================================================================================
225!  subroutine Agrif_detect_All
226!
227!> Detects areas to be refined.
228!---------------------------------------------------------------------------------------------------
229recursive subroutine Agrif_detect_all ( g )
230!---------------------------------------------------------------------------------------------------
231    TYPE(Agrif_Grid),  pointer  :: g        !< Pointer on the current grid
232!
233    Type(Agrif_PGrid), pointer  :: parcours ! Pointer for the recursive procedure
234    integer, DIMENSION(3)       :: size
235    integer                     :: i
236    real                        :: g_eps
237!
238    parcours => g % child_list % first
239!
240!   To be positioned on the finer grids of the grid hierarchy
241!
242    do while (associated(parcours))
243        call Agrif_detect_all(parcours % gr)
244        parcours => parcours % next
245    enddo
246!
247    g_eps = huge(1.)
248    do i = 1,Agrif_Probdim
249        g_eps = min(g_eps, g % Agrif_dx(i))
250    enddo
251!
252    g_eps = g_eps / 100.
253!
254    if ( Agrif_Probdim == 1 ) g%tabpoint1D = 0
255    if ( Agrif_Probdim == 2 ) g%tabpoint2D = 0
256    if ( Agrif_Probdim == 3 ) g%tabpoint3D = 0
257!
258    do i = 1,Agrif_Probdim
259        if ( g%Agrif_dx(i)/Agrif_coeffref(i) < (Agrif_mind(i)-g_eps) ) return
260    enddo
261!
262    call Agrif_instance(g)
263!
264!   Detection (Agrif_detect is a users routine)
265!
266    do i = 1,Agrif_Probdim
267        size(i) = g % nb(i) + 1
268    enddo
269!
270    SELECT CASE (Agrif_Probdim)
271    CASE (1)
272        call Agrif_detect(g%tabpoint1D,size)
273    CASE (2)
274        call Agrif_detect(g%tabpoint2D,size)
275    CASE (3)
276        call Agrif_detect(g%tabpoint3D,size)
277    END SELECT
278!
279!   Addition of the areas detected on the child grids
280!
281    parcours => g % child_list % first
282!
283    do while (associated(parcours))
284        call Agrif_Add_detected_areas(g,parcours % gr)
285        parcours => parcours % next
286    enddo
287!---------------------------------------------------------------------------------------------------
288end subroutine Agrif_detect_all
289!===================================================================================================
290!
291!===================================================================================================
292!  subroutine Agrif_Add_detected_areas
293!
294!> Adds on the parent grid the areas detected on its child grids
295!---------------------------------------------------------------------------------------------------
296subroutine Agrif_Add_detected_areas ( parentgrid, childgrid )
297!---------------------------------------------------------------------------------------------------
298    Type(Agrif_Grid), pointer   :: parentgrid
299    Type(Agrif_Grid), pointer   :: childgrid
300!
301    integer :: i,j,k
302!
303    do i = 1,childgrid%nb(1)+1
304        if ( Agrif_Probdim == 1 ) then
305            if (childgrid%tabpoint1D(i)==1) then
306                parentgrid%tabpoint1D(childgrid%ix(1)+(i-1)/Agrif_Coeffref(1)) = 1
307            endif
308        else
309            do j=1,childgrid%nb(2)+1
310                if (Agrif_Probdim==2) then
311                    if (childgrid%tabpoint2D(i,j)==1) then
312                        parentgrid%tabpoint2D(                       &
313                            childgrid%ix(1)+(i-1)/Agrif_Coeffref(1), &
314                            childgrid%ix(2)+(j-1)/Agrif_Coeffref(2)) = 1
315                    endif
316                else
317                    do k=1,childgrid%nb(3)+1
318                        if (childgrid%tabpoint3D(i,j,k)==1) then
319                            parentgrid%tabpoint3D(                       &
320                                childgrid%ix(1)+(i-1)/Agrif_Coeffref(1), &
321                                childgrid%ix(2)+(j-1)/Agrif_Coeffref(2), &
322                                childgrid%ix(3)+(k-1)/Agrif_Coeffref(3)) = 1
323                        endif
324                    enddo
325                endif
326            enddo
327        endif
328    enddo
329!---------------------------------------------------------------------------------------------------
330end subroutine Agrif_Add_detected_areas
331!===================================================================================================
332!
333!===================================================================================================
334!  subroutine Agrif_Free_before_All
335!---------------------------------------------------------------------------------------------------
336recursive subroutine Agrif_Free_before_All ( gridlist )
337!---------------------------------------------------------------------------------------------------
338    Type(Agrif_Grid_List), intent(inout)    :: gridlist !< Grid list
339!
340    Type(Agrif_PGrid), pointer   :: parcours ! Pointer for the recursive procedure
341!
342    parcours => gridlist % first
343!
344    do while (associated(parcours))
345!
346        if (.not. parcours%gr%fixed) then
347            call Agrif_Free_data_before(parcours%gr)
348            parcours % gr % oldgrid = .TRUE.
349        endif
350!
351        call Agrif_Free_before_all (parcours % gr % child_list)
352!
353        parcours => parcours % next
354!
355    enddo
356!---------------------------------------------------------------------------------------------------
357end subroutine Agrif_Free_before_All
358!===================================================================================================
359!
360!===================================================================================================
361!  subroutine Agrif_Save_All
362!---------------------------------------------------------------------------------------------------
363recursive subroutine Agrif_Save_All ( gridlist )
364!---------------------------------------------------------------------------------------------------
365    type(Agrif_Grid_List), intent(inout)    :: gridlist !< Grid list
366!
367    type(Agrif_PGrid), pointer   :: parcours ! Pointer for the recursive procedure
368!
369    parcours => gridlist % first
370!
371    do while (associated(parcours))
372!
373        if (.not. parcours%gr%fixed) then
374            call Agrif_Instance(parcours%gr)
375            call Agrif_Before_Regridding()
376            parcours % gr % oldgrid = .TRUE.
377        endif
378!
379        call Agrif_Save_All(parcours % gr % child_list)
380!
381        parcours => parcours % next
382!
383      enddo
384!---------------------------------------------------------------------------------------------------
385end subroutine Agrif_Save_All
386!===================================================================================================
387!
388!===================================================================================================
389!  subroutine Agrif_Free_after_All
390!---------------------------------------------------------------------------------------------------
391recursive subroutine Agrif_Free_after_All ( gridlist )
392!---------------------------------------------------------------------------------------------------
393    Type(Agrif_Grid_List), intent(inout)    :: gridlist       !< Grid list to free
394!
395    Type(Agrif_PGrid), pointer  :: parcours ! Pointer for the recursive proced
396    Type(Agrif_PGrid), pointer  :: preparcours
397    Type(Agrif_PGrid), pointer  :: preparcoursini
398!
399    allocate(preparcours)
400!
401    preparcoursini => preparcours
402!
403    nullify(preparcours % gr)
404!
405    preparcours % next => gridlist % first
406    parcours => gridlist % first
407!
408    do while (associated(parcours))
409!
410        if ( (.NOT. parcours%gr % fixed) .AND. (parcours%gr % oldgrid) ) then
411            call Agrif_Free_data_after(parcours%gr)
412        endif
413!
414        call Agrif_Free_after_all( parcours%gr % child_list )
415!
416        if (parcours % gr % oldgrid) then
417            deallocate(parcours % gr)
418            preparcours % next => parcours % next
419            deallocate(parcours)
420            parcours => preparcours % next
421        else
422            preparcours => preparcours % next
423            parcours => parcours % next
424        endif
425!
426    enddo
427!
428    deallocate(preparcoursini)
429!---------------------------------------------------------------------------------------------------
430end subroutine Agrif_Free_after_All
431!===================================================================================================
432!
433!===================================================================================================
434!  subroutine Agrif_Integrate
435!
436!> Manages the time integration of the grid hierarchy.
437!! Recursive subroutine and call on subroutines Agrif_Init::Agrif_Instance and #Agrif_Step
438!---------------------------------------------------------------------------------------------------
439recursive subroutine Agrif_Integrate ( g, procname )
440!---------------------------------------------------------------------------------------------------
441    type(Agrif_Grid), pointer   :: g        !< Pointer on the current grid
442    procedure(step_proc)        :: procname !< Subroutine to call on each grid
443!
444    type(Agrif_PGrid), pointer  :: parcours ! Pointer for the recursive procedure
445    integer                     :: nbt      ! Number of time steps of the current grid
446    integer                     :: i, k
447!
448!   Instanciation of the variables of the current grid
449!    if ( g % fixedrank /= 0 ) then
450        call Agrif_Instance(g)
451!    endif
452!
453!   One step on the current grid
454!
455    call procname ()
456!
457!   Number of time steps on the current grid
458!
459    g%ngridstep = g % ngridstep + 1
460    parcours => g % child_list % first
461!
462!   Recursive procedure for the time integration of the grid hierarchy
463    do while (associated(parcours))
464!
465!       Instanciation of the variables of the current grid
466        call Agrif_Instance(parcours % gr)
467!
468!       Number of time steps
469        nbt = 1
470        do i = 1,Agrif_Probdim
471            nbt = max(nbt, parcours % gr % timeref(i))
472        enddo
473!
474        do k = 1,nbt
475            call Agrif_Integrate(parcours % gr, procname)
476        enddo
477!
478        parcours => parcours % next
479!
480    enddo
481!---------------------------------------------------------------------------------------------------
482end subroutine Agrif_Integrate
483!===================================================================================================
484!
485!===================================================================================================
486!  subroutine Agrif_Integrate_Parallel
487!
488!> Manages the time integration of the grid hierarchy in parallel
489!! Recursive subroutine and call on subroutines Agrif_Init::Agrif_Instance and #Agrif_Step
490!---------------------------------------------------------------------------------------------------
491recursive subroutine Agrif_Integrate_Parallel ( g, procname )
492!---------------------------------------------------------------------------------------------------
493    type(Agrif_Grid), pointer   :: g        !< Pointer on the current grid
494    procedure(step_proc)        :: procname !< Subroutine to call on each grid
495!
496#if defined AGRIF_MPI
497    type(Agrif_PGrid), pointer  :: gridp    ! Pointer for the recursive procedure
498    integer                     :: nbt      ! Number of time steps of the current grid
499    integer                     :: i, k, is
500!
501!   Instanciation of the variables of the current grid
502    if ( g % fixedrank /= 0 ) then
503        call Agrif_Instance(g)
504    endif
505!
506! One step on the current grid
507    call procname ()
508!
509! Number of time steps on the current grid
510    g % ngridstep = g % ngridstep + 1
511!
512! Continue only if the grid has defined sequences of child integrations.
513    if ( .not. associated(g % child_seq) ) return
514!
515    do is = 1, g % child_seq % nb_seqs
516!
517!     For each sequence, a given processor does integrate only on grid.
518        gridp => Agrif_seq_select_child(g,is)
519!
520!     Instanciation of the variables of the current grid
521        call Agrif_Instance(gridp % gr)
522!
523!     Number of time steps
524        nbt = 1
525        do i = 1,Agrif_Probdim
526            nbt = max(nbt, gridp % gr % timeref(i))
527        enddo
528!
529        do k = 1,nbt
530            call Agrif_Integrate_Parallel(gridp % gr, procname)
531        enddo
532!
533    enddo
534#else
535    call Agrif_Integrate( g, procname )
536#endif
537!---------------------------------------------------------------------------------------------------
538end subroutine Agrif_Integrate_Parallel
539!===================================================================================================
540!
541!
542!===================================================================================================
543!  subroutine Agrif_Integrate_ChildGrids
544!
545!> Manages the time integration of the grid hierarchy.
546!! Call the subroutine procname on each child grid of the current grid
547!---------------------------------------------------------------------------------------------------
548recursive subroutine Agrif_Integrate_ChildGrids ( procname )
549!---------------------------------------------------------------------------------------------------
550    procedure(step_proc)        :: procname !< Subroutine to call on each grid
551!
552    type(Agrif_PGrid), pointer  :: parcours ! Pointer for the recursive procedure
553    integer                     :: nbt      ! Number of time steps of the current grid
554    integer                     :: i, k, is
555    type(Agrif_Grid)                    , pointer :: save_grid
556    type(Agrif_PGrid), pointer  :: gridp    ! Pointer for the recursive procedure
557   
558    save_grid => Agrif_Curgrid
559
560! Number of time steps on the current grid
561    save_grid % ngridstep = save_grid % ngridstep + 1
562   
563#ifdef AGRIF_MPI
564    if ( .not. Agrif_Parallel_sisters ) then
565#endif
566    parcours => save_grid % child_list % first
567!
568!   Recursive procedure for the time integration of the grid hierarchy
569    do while (associated(parcours))
570!
571!       Instanciation of the variables of the current grid
572        call Agrif_Instance(parcours % gr)
573!
574!       Number of time steps
575        nbt = 1
576        do i = 1,Agrif_Probdim
577            nbt = max(nbt, parcours % gr % timeref(i))
578        enddo
579!
580        do k = 1,nbt
581            call procname()
582        enddo
583!
584        parcours => parcours % next
585!
586    enddo
587
588#ifdef AGRIF_MPI
589    else
590! Continue only if the grid has defined sequences of child integrations.
591    if ( .not. associated(save_grid % child_seq) ) return
592!
593    do is = 1, save_grid % child_seq % nb_seqs
594!
595!     For each sequence, a given processor does integrate only on grid.
596        gridp => Agrif_seq_select_child(save_grid,is)
597!
598!     Instanciation of the variables of the current grid
599        call Agrif_Instance(gridp % gr)
600!
601!     Number of time steps
602        nbt = 1
603        do i = 1,Agrif_Probdim
604            nbt = max(nbt, gridp % gr % timeref(i))
605        enddo
606!
607        do k = 1,nbt
608            call procname()
609        enddo
610!
611    enddo
612    endif
613#endif
614
615    call Agrif_Instance(save_grid)
616   
617!---------------------------------------------------------------------------------------------------
618end subroutine Agrif_Integrate_ChildGrids
619!===================================================================================================
620!===================================================================================================
621!  subroutine Agrif_Integrate_Child
622!
623!> Manages the time integration of the grid hierarchy.
624!! Recursive subroutine and call on subroutines Agrif_Instance & Agrif_Step.
625!---------------------------------------------------------------------------------------------------
626recursive subroutine Agrif_Integrate_Child ( g, procname )
627!---------------------------------------------------------------------------------------------------
628    type(Agrif_Grid), pointer   :: g        !< Pointer on the current grid
629    procedure(step_proc)        :: procname !< Subroutine to call on each grid
630!
631    type(Agrif_PGrid), pointer  :: parcours ! Pointer for the recursive procedure
632!
633!   One step on the current grid
634!
635    call procname ()
636!
637!   Number of time steps on the current grid
638!
639    parcours => g % child_list % first
640!
641!   Recursive procedure for the time integration of the grid hierarchy
642    do while (associated(parcours))
643!
644!       Instanciation of the variables of the current grid
645        call Agrif_Instance(parcours % gr)
646        call Agrif_Integrate_Child (parcours % gr, procname)
647        parcours => parcours % next
648!
649    enddo
650!---------------------------------------------------------------------------------------------------
651end subroutine Agrif_Integrate_Child
652!===================================================================================================
653!
654!===================================================================================================
655!  subroutine Agrif_Integrate_Child_Parallel
656!
657!> Manages the time integration of the grid hierarchy.
658!! Recursive subroutine and call on subroutines Agrif_Instance & Agrif_Step.
659!---------------------------------------------------------------------------------------------------
660recursive subroutine Agrif_Integrate_Child_Parallel ( g, procname )
661!---------------------------------------------------------------------------------------------------
662    type(Agrif_Grid), pointer   :: g        !< Pointer on the current grid
663    procedure(step_proc)        :: procname !< Subroutine to call on each grid
664!
665#if defined AGRIF_MPI
666    type(Agrif_PGrid), pointer  :: gridp    ! Pointer for the recursive procedure
667    integer                     :: is
668!
669! Instanciation of the variables of the current grid
670    call Agrif_Instance(g)
671!
672! One step on the current grid
673    call procname ()
674!
675! Continue only if the grid has defined sequences of child integrations.
676    if ( .not. associated(g % child_seq) ) return
677!
678    do is = 1, g % child_seq % nb_seqs
679!
680!     For each sequence, a given processor does integrate only on grid.
681        gridp => Agrif_seq_select_child(g,is)
682        call Agrif_Integrate_Child_Parallel(gridp % gr, procname)
683!
684    enddo
685!
686    call Agrif_Instance(g)
687#else
688    call Agrif_Integrate_Child( g, procname )
689#endif
690!---------------------------------------------------------------------------------------------------
691end subroutine Agrif_Integrate_Child_Parallel
692!===================================================================================================
693!
694!===================================================================================================
695!  subroutine Agrif_Init_Grids
696!
697!> Initializes the root coarse grid pointed by Agrif_Mygrid. It is called in the main program.
698!---------------------------------------------------------------------------------------------------
699subroutine Agrif_Init_Grids ( procname1, procname2 )
700!---------------------------------------------------------------------------------------------------
701    procedure(typedef_proc), optional   :: procname1 !< (Default: Agrif_probdim_modtype_def)
702    procedure(alloc_proc),   optional   :: procname2 !< (Default: Agrif_Allocationcalls)
703!
704    integer :: i, ierr_allocate, nunit
705    integer :: is_coarse, rhox,rhoy,rhoz,rhot
706    logical :: BEXIST
707!
708    if (present(procname1)) Then
709        call procname1()
710    else
711        call Agrif_probdim_modtype_def()
712    endif
713!
714
715! TEST FOR COARSE GRID (GRAND MOTHER GRID) in AGRIF_FixedGrids.in
716    nunit = Agrif_Get_Unit()
717    open(nunit, file='AGRIF_FixedGrids.in', form='formatted', status="old", ERR=98)
718    if (Agrif_Probdim == 3) then
719       read(nunit,*) is_coarse, rhox, rhoy, rhoz, rhot
720    elseif (Agrif_Probdim == 2) then
721       read(nunit,*) is_coarse, rhox, rhoy, rhot
722    elseif (Agrif_Probdim == 2) then
723       read(nunit,*) is_coarse, rhox, rhot
724    endif
725    if (is_coarse == -1) then
726       agrif_coarse = .TRUE.
727       if (Agrif_Probdim == 3) then
728          coarse_spaceref(1:3)=(/rhox,rhoy,rhoz/)
729       elseif (Agrif_Probdim == 2) then
730          coarse_spaceref(1:2)=(/rhox,rhoy/)
731       elseif (Agrif_Probdim == 2) then
732          coarse_spaceref(1:1)=(/rhox/)
733       endif
734       coarse_timeref(1:Agrif_Probdim) = rhot
735    endif
736    close(nunit)
737   
738    Agrif_UseSpecialValue = .FALSE.
739    Agrif_UseSpecialValueFineGrid = .FALSE.
740    Agrif_SpecialValue = 0.
741    Agrif_SpecialValueFineGrid = 0.
742!
743    allocate(Agrif_Mygrid)
744    allocate(Agrif_OldMygrid)
745!
746!   Space and time refinement factors are set to 1 on the root grid
747!
748    do i = 1,Agrif_Probdim
749        Agrif_Mygrid % spaceref(i) = coarse_spaceref(i)
750        Agrif_Mygrid % timeref(i)  = coarse_timeref(i)
751    enddo
752!
753!   Initialization of the number of time steps
754    Agrif_Mygrid % ngridstep = 0
755    Agrif_Mygrid % grid_id   = 0
756!
757!   No parent grid for the root coarse grid
758    nullify(Agrif_Mygrid % parent)
759!
760!   Initialization of the minimum positions, global abscissa and space steps
761    do i = 1, Agrif_Probdim
762        Agrif_Mygrid % ix(i) = 1
763        Agrif_Mygrid % Agrif_x(i)  = 0.
764        Agrif_Mygrid % Agrif_dx(i) = 1./Agrif_Mygrid % spaceref(i)
765        Agrif_Mygrid % Agrif_dt(i) = 1./Agrif_Mygrid % timeref(i)
766!       Borders of the root coarse grid
767        Agrif_Mygrid % NearRootBorder(i) = .true.
768        Agrif_Mygrid % DistantRootBorder(i) = .true.
769    enddo
770!
771!   The root coarse grid is a fixed grid
772    Agrif_Mygrid % fixed = .TRUE.
773!   Level of the root grid
774    Agrif_Mygrid % level = 0
775!   Maximum level in the hierarchy
776    Agrif_MaxLevelLoc = 0
777!
778!   Number of the grid pointed by Agrif_Mygrid (root coarse grid)
779    Agrif_Mygrid % rank = 1
780!
781!   Number of the root grid as a fixed grid
782    Agrif_Mygrid % fixedrank = 0
783!
784!   Initialization of some fields of the root grid variables
785    ierr_allocate = 0
786    if( Agrif_NbVariables(0) > 0 ) allocate(Agrif_Mygrid % tabvars(Agrif_NbVariables(0)),stat=ierr_allocate)
787    if( Agrif_NbVariables(1) > 0 ) allocate(Agrif_Mygrid % tabvars_c(Agrif_NbVariables(1)),stat=ierr_allocate)
788    if( Agrif_NbVariables(2) > 0 ) allocate(Agrif_Mygrid % tabvars_r(Agrif_NbVariables(2)),stat=ierr_allocate)
789    if( Agrif_NbVariables(3) > 0 ) allocate(Agrif_Mygrid % tabvars_l(Agrif_NbVariables(3)),stat=ierr_allocate)
790    if( Agrif_NbVariables(4) > 0 ) allocate(Agrif_Mygrid % tabvars_i(Agrif_NbVariables(4)),stat=ierr_allocate)
791    if (ierr_allocate /= 0) THEN
792      STOP "*** ERROR WHEN ALLOCATING TABVARS ***"
793    endif
794!
795!   Initialization of the other fields of the root grid variables (number of
796!   cells, positions, number and type of its dimensions, ...)
797    call Agrif_Instance(Agrif_Mygrid)
798    call Agrif_Set_numberofcells(Agrif_Mygrid)
799!
800!   Allocation of the array containing the values of the grid variables
801    call Agrif_Allocation(Agrif_Mygrid, procname2)
802    call Agrif_initialisations(Agrif_Mygrid)
803!
804!   Total number of fixed grids
805    Agrif_nbfixedgrids = 0
806   
807! If a grand mother grid is declared
808
809    if (agrif_coarse) then
810      allocate(Agrif_Coarsegrid)
811
812      Agrif_Coarsegrid % ngridstep = 0
813      Agrif_Coarsegrid % grid_id   = -9999
814     
815    do i = 1, Agrif_Probdim
816        Agrif_Coarsegrid%spaceref(i) = coarse_spaceref(i)
817        Agrif_Coarsegrid%timeref(i) = coarse_timeref(i)
818        Agrif_Coarsegrid % ix(i) = 1
819        Agrif_Coarsegrid % Agrif_x(i)  = 0.
820        Agrif_Coarsegrid % Agrif_dx(i) = 1.
821        Agrif_Coarsegrid % Agrif_dt(i) = 1.
822!       Borders of the root coarse grid
823        Agrif_Coarsegrid % NearRootBorder(i) = .true.
824        Agrif_Coarsegrid % DistantRootBorder(i) = .true.
825        Agrif_Coarsegrid % nb(i) =Agrif_mygrid%nb(i) / coarse_spaceref(i)
826    enddo     
827
828!   The root coarse grid is a fixed grid
829    Agrif_Coarsegrid % fixed = .TRUE.
830!   Level of the root grid
831    Agrif_Coarsegrid % level = -1
832   
833    Agrif_Coarsegrid % grand_mother_grid = .true.
834
835!   Number of the grid pointed by Agrif_Mygrid (root coarse grid)
836    Agrif_Coarsegrid % rank = -9999
837!
838!   Number of the root grid as a fixed grid
839    Agrif_Coarsegrid % fixedrank = -9999
840   
841      Agrif_Mygrid%parent => Agrif_Coarsegrid
842     
843! Not used but required to prevent seg fault
844      Agrif_Coarsegrid%parent => Agrif_Mygrid
845     
846      call Agrif_Create_Var(Agrif_Coarsegrid)
847
848! Reset to null
849      Nullify(Agrif_Coarsegrid%parent)
850     
851      Agrif_Coarsegrid%child_list%nitems = 1
852      allocate(Agrif_Coarsegrid%child_list%first)
853      allocate(Agrif_Coarsegrid%child_list%last)
854      Agrif_Coarsegrid%child_list%first%gr => Agrif_Mygrid
855      Agrif_Coarsegrid%child_list%last%gr => Agrif_Mygrid
856
857    endif
858   
859    return
860
86198  INQUIRE(FILE='AGRIF_FixedGrids.in',EXIST=BEXIST)
862    if (.not. BEXIST) then
863        print*,'ERROR : File AGRIF_FixedGrids.in not found.'
864        STOP
865    else
866        print*,'Error opening file AGRIF_FixedGrids.in'
867        STOP
868    endif
869   
870!---------------------------------------------------------------------------------------------------
871end subroutine Agrif_Init_Grids
872!===================================================================================================
873!
874!===================================================================================================
875!  subroutine Agrif_Deallocation
876!
877!> Deallocates all data arrays.
878!---------------------------------------------------------------------------------------------------
879subroutine Agrif_Deallocation
880!---------------------------------------------------------------------------------------------------
881    integer                         :: nb
882    type(Agrif_Variable), pointer   :: var
883    type(Agrif_Variable_c), pointer   :: var_c
884    type(Agrif_Variable_l), pointer   :: var_l
885    type(Agrif_Variable_i), pointer   :: var_i
886!
887    do nb = 1,Agrif_NbVariables(0)
888!
889        var => Agrif_Mygrid % tabvars(nb)
890!
891        if ( allocated(var % array1) ) deallocate(var % array1)
892        if ( allocated(var % array2) ) deallocate(var % array2)
893        if ( allocated(var % array3) ) deallocate(var % array3)
894        if ( allocated(var % array4) ) deallocate(var % array4)
895        if ( allocated(var % array5) ) deallocate(var % array5)
896        if ( allocated(var % array6) ) deallocate(var % array6)
897!
898        if ( allocated(var % sarray1) ) deallocate(var % sarray1)
899        if ( allocated(var % sarray2) ) deallocate(var % sarray2)
900        if ( allocated(var % sarray3) ) deallocate(var % sarray3)
901        if ( allocated(var % sarray4) ) deallocate(var % sarray4)
902        if ( allocated(var % sarray5) ) deallocate(var % sarray5)
903        if ( allocated(var % sarray6) ) deallocate(var % sarray6)
904!
905        if ( allocated(var % darray1) ) deallocate(var % darray1)
906        if ( allocated(var % darray2) ) deallocate(var % darray2)
907        if ( allocated(var % darray3) ) deallocate(var % darray3)
908        if ( allocated(var % darray4) ) deallocate(var % darray4)
909        if ( allocated(var % darray5) ) deallocate(var % darray5)
910        if ( allocated(var % darray6) ) deallocate(var % darray6)
911!
912    enddo
913!
914    do nb = 1,Agrif_NbVariables(1)
915!
916        var_c => Agrif_Mygrid % tabvars_c(nb)
917!
918        if ( allocated(var_c % carray1) ) deallocate(var_c % carray1)
919        if ( allocated(var_c % carray2) ) deallocate(var_c % carray2)
920!
921    enddo
922
923    do nb = 1,Agrif_NbVariables(3)
924!
925        var_l => Agrif_Mygrid % tabvars_l(nb)
926!
927        if ( allocated(var_l % larray1) ) deallocate(var_l % larray1)
928        if ( allocated(var_l % larray2) ) deallocate(var_l % larray2)
929        if ( allocated(var_l % larray3) ) deallocate(var_l % larray3)
930        if ( allocated(var_l % larray4) ) deallocate(var_l % larray4)
931        if ( allocated(var_l % larray5) ) deallocate(var_l % larray5)
932        if ( allocated(var_l % larray6) ) deallocate(var_l % larray6)
933!
934    enddo
935!
936    do nb = 1,Agrif_NbVariables(4)
937!
938        var_i => Agrif_Mygrid % tabvars_i(nb)
939!
940        if ( allocated(var_i % iarray1) ) deallocate(var_i % iarray1)
941        if ( allocated(var_i % iarray2) ) deallocate(var_i % iarray2)
942        if ( allocated(var_i % iarray3) ) deallocate(var_i % iarray3)
943        if ( allocated(var_i % iarray4) ) deallocate(var_i % iarray4)
944        if ( allocated(var_i % iarray5) ) deallocate(var_i % iarray5)
945        if ( allocated(var_i % iarray6) ) deallocate(var_i % iarray6)
946!
947    enddo
948!
949    if ( allocated(Agrif_Mygrid % tabvars) )   deallocate(Agrif_Mygrid % tabvars)
950    if ( allocated(Agrif_Mygrid % tabvars_c) ) deallocate(Agrif_Mygrid % tabvars_c)
951    if ( allocated(Agrif_Mygrid % tabvars_r) ) deallocate(Agrif_Mygrid % tabvars_r)
952    if ( allocated(Agrif_Mygrid % tabvars_l) ) deallocate(Agrif_Mygrid % tabvars_l)
953    if ( allocated(Agrif_Mygrid % tabvars_i) ) deallocate(Agrif_Mygrid % tabvars_i)
954    deallocate(Agrif_Mygrid)
955!---------------------------------------------------------------------------------------------------
956end subroutine Agrif_Deallocation
957!===================================================================================================
958!
959!===================================================================================================
960!  subroutine Agrif_Step_adj
961!
962!> creates the grid hierarchy and manages the backward time integration procedure.
963!> It is called in the main program.
964!> calls subroutines #Agrif_Regrid and #Agrif_Integrate_adj.
965!---------------------------------------------------------------------------------------------------
966subroutine Agrif_Step_adj ( procname )
967!---------------------------------------------------------------------------------------------------
968    procedure(step_proc)    :: procname     !< Subroutine to call on each grid
969!
970!   Creation and initialization of the grid hierarchy
971!
972!   Set the clustering variables
973    call Agrif_clustering_def()
974!
975    if ( Agrif_USE_ONLY_FIXED_GRIDS .EQ. 1 ) then
976!
977        if (Agrif_Mygrid % ngridstep == 0) then
978            if (.not.Agrif_regrid_has_been_done ) then
979                call Agrif_Regrid()
980            endif
981            call Agrif_Instance(Agrif_Mygrid)
982        endif
983!
984    else
985!
986        if (mod(Agrif_Mygrid % ngridstep, Agrif_Regridding) == 0) then
987            call Agrif_Regrid()
988            call Agrif_Instance(Agrif_Mygrid)
989        endif
990!
991    endif
992!
993!   Time integration of the grid hierarchy
994!
995    call Agrif_Integrate_adj (Agrif_Mygrid,procname)
996!
997    if ( Agrif_Mygrid % child_list % nitems > 0 ) call Agrif_Instance(Agrif_Mygrid)
998!
999!---------------------------------------------------------------------------------------------------
1000end subroutine Agrif_Step_adj
1001!===================================================================================================
1002!
1003!===================================================================================================
1004!  subroutine Agrif_Integrate_adj
1005!
1006!> Manages the backward time integration of the grid hierarchy.
1007!! Recursive subroutine and call on subroutines Agrif_Init::Agrif_Instance and #Agrif_Step_adj
1008!---------------------------------------------------------------------------------------------------
1009recursive subroutine Agrif_Integrate_adj ( g, procname )
1010!---------------------------------------------------------------------------------------------------
1011    type(Agrif_Grid), pointer   :: g        !< Pointer on the current grid
1012    procedure(step_proc)        :: procname !< Subroutine to call on each grid
1013!
1014    type(Agrif_pgrid), pointer :: parcours ! pointer for the recursive procedure
1015    integer                    :: nbt      ! Number of time steps of the current grid
1016    integer                    :: k
1017!
1018!   Instanciation of the variables of the current grid
1019    if ( g%fixedrank /= 0 ) then
1020        call Agrif_Instance(g)
1021    endif
1022!
1023!   Number of time steps on the current grid
1024!
1025    g%ngridstep = g % ngridstep + 1
1026    parcours => g % child_list % first
1027!
1028!   Recursive procedure for the time integration of the grid hierarchy
1029    do while (associated(parcours))
1030!
1031!   Instanciation of the variables of the current grid
1032        call Agrif_Instance(parcours % gr)
1033!
1034!       Number of time steps
1035        nbt = 1
1036        do k = 1,Agrif_Probdim
1037            nbt = max(nbt, parcours % gr % timeref(k))
1038        enddo
1039!
1040        do k = nbt,1,-1
1041            call Agrif_Integrate_adj(parcours % gr, procname)
1042        enddo
1043!
1044        parcours => parcours % next
1045!
1046    enddo
1047!
1048    if ( g % child_list % nitems > 0 )  call Agrif_Instance(g)
1049!
1050!   One step on the current grid
1051    call procname ()
1052!
1053end subroutine Agrif_Integrate_adj
1054!===================================================================================================
1055!
1056!===================================================================================================
1057!  subroutine Agrif_Step_Child_adj
1058!
1059!> Apply 'procname' to each grid of the hierarchy from Child to Parent
1060!---------------------------------------------------------------------------------------------------
1061subroutine Agrif_Step_Child_adj ( procname )
1062!---------------------------------------------------------------------------------------------------
1063    procedure(step_proc)        :: procname !< Subroutine to call on each grid
1064!
1065    call Agrif_Integrate_Child_adj(Agrif_Mygrid,procname)
1066!
1067    if ( Agrif_Mygrid % child_list % nitems > 0 ) call Agrif_Instance(Agrif_Mygrid)
1068!
1069end subroutine Agrif_Step_Child_adj
1070!===================================================================================================
1071!
1072!===================================================================================================
1073!  subroutine Agrif_Integrate_Child_adj
1074!
1075!> Manages the backward time integration of the grid hierarchy.
1076!! Recursive subroutine and call on subroutines Agrif_Init::Agrif_Instance & Agrif_Step_adj.
1077!---------------------------------------------------------------------------------------------------
1078recursive subroutine Agrif_Integrate_Child_adj ( g, procname )
1079!---------------------------------------------------------------------------------------------------
1080    type(Agrif_Grid),pointer   :: g          !< Pointer on the current grid
1081    procedure(step_proc)       :: procname   !< Subroutine to call on each grid
1082!
1083    type(Agrif_PGrid),pointer  :: parcours   !< Pointer for the recursive procedure
1084!
1085    parcours => g % child_list % first
1086!
1087!   Recursive procedure for the time integration of the grid hierarchy
1088    do while (associated(parcours))
1089!
1090!       Instanciation of the variables of the current grid
1091        call Agrif_Instance(parcours % gr)
1092        call Agrif_Integrate_Child_adj(parcours % gr, procname)
1093!
1094       parcours => parcours % next
1095!
1096    enddo
1097    if ( g % child_list % nitems > 0 )  call Agrif_Instance(g)
1098!
1099!   One step on the current grid
1100    call procname()
1101!---------------------------------------------------------------------------------------------------
1102end subroutine Agrif_Integrate_Child_adj
1103!===================================================================================================
1104!
1105!===================================================================================================
1106
1107end module Agrif_Util
Note: See TracBrowser for help on using the repository browser.