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

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

Last change on this file was 14107, checked in by nicolasmartin, 3 years ago

Reintegration of dev_r12970_AGRIF_CMEMS to AGRIF/dev

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