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

source: vendors/AGRIF/CMEMS_2020/AGRIF_FILES/modseq.F90 @ 10087

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

update AGRIF library

  • Property svn:keywords set to Id
File size: 24.7 KB
Line 
1module Agrif_seq
2!
3    use Agrif_Init
4    use Agrif_Procs
5    use Agrif_Grids
6    !use Agrif_Arrays
7!
8    implicit none
9!
10contains
11!
12#if defined AGRIF_MPI
13!===================================================================================================
14function Agrif_seq_allocate_list ( nb_seqs ) result( seqlist )
15!---------------------------------------------------------------------------------------------------
16    integer, intent(in)     :: nb_seqs
17!
18    type(Agrif_Sequence_List), pointer :: seqlist
19!
20    allocate(seqlist)
21    seqlist % nb_seqs = nb_seqs
22    allocate(seqlist % sequences(1:nb_seqs))
23!---------------------------------------------------------------------------------------------------
24end function Agrif_seq_allocate_list
25!===================================================================================================
26!
27!===================================================================================================
28subroutine Agrif_seq_add_grid ( seqlist, seq_num, grid )
29!---------------------------------------------------------------------------------------------------
30    type(Agrif_Sequence_List), intent(inout)    :: seqlist
31    integer,                   intent(in)       :: seq_num
32    type(Agrif_Grid), pointer, intent(in)       :: grid
33!
34    call Agrif_gl_append(seqlist % sequences(seq_num) % gridlist, grid )
35!---------------------------------------------------------------------------------------------------
36end subroutine Agrif_seq_add_grid
37!===================================================================================================
38!
39!===================================================================================================
40recursive subroutine Agrif_seq_init_sequences ( grid )
41!---------------------------------------------------------------------------------------------------
42    type(Agrif_Grid), intent(inout)  :: grid
43!
44    type(Agrif_PGrid), pointer  :: gp
45!
46#if defined AGRIF_MPI
47!
48! Build list of required procs for each child grid
49    gp => grid % child_list % first
50    do while ( associated( gp ) )
51        call Agrif_seq_build_required_proclist( gp % gr )
52        gp => gp % next
53    enddo
54!
55! Create integration sequences for the current grid
56    call Agrif_seq_create_proc_sequences( grid )
57    call Agrif_seq_allocate_procs_to_childs( grid )
58!
59! Create new communicators for sequences
60    call Agrif_seq_create_communicators( grid )
61!
62#endif
63!---------------------------------------------------------------------------------------------------
64end subroutine Agrif_seq_init_sequences
65!===================================================================================================
66!
67!===================================================================================================
68subroutine Agrif_seq_build_required_proclist ( grid )
69!---------------------------------------------------------------------------------------------------
70    type(Agrif_Grid), intent(inout) :: grid
71!
72    type(Agrif_Grid),      pointer  :: parent_grid
73    type(Agrif_Rectangle), pointer  :: grid_rect
74    type(Agrif_Proc_p),    pointer  :: proc_rect
75    type(Agrif_Proc),      pointer  :: proc
76    logical     :: proc_is_required
77    integer     :: i
78!
79    if ( grid % fixedrank == 0 ) then
80!       grid is the Root
81        if ( grid % required_proc_list % nitems == 0 ) then
82            print*, "### Error Agrif_seq_build_required_proclist: empty proc list."
83            print*, "# -> You should check if Agrif_Init_ProcList() is actually working."
84            stop
85        endif
86        return
87    endif
88!
89    parent_grid => grid % parent
90    grid_rect   => grid % rect_in_parent
91    proc_rect   => parent_grid % proc_def_list % first
92!
93
94    do while ( associated( proc_rect ) )
95
96!
97        proc => proc_rect % proc
98!
99        proc_is_required = .true.
100        do i = 1,Agrif_Probdim
101            proc_is_required = ( proc_is_required             .and. &
102                    ( grid_rect % imin(i) <= proc % imax(i) ) .and. &
103                    ( grid_rect % imax(i) >= proc % imin(i) ) )
104        enddo
105!
106        if ( proc_is_required ) then
107            call Agrif_pl_append(grid % required_proc_list, proc)
108        endif
109        proc_rect => proc_rect % next
110!
111    enddo
112!---------------------------------------------------------------------------------------------------
113end subroutine Agrif_seq_build_required_proclist
114!===================================================================================================
115!
116!===================================================================================================
117subroutine Agrif_seq_create_proc_sequences ( grid )
118!---------------------------------------------------------------------------------------------------
119    type(Agrif_Grid), intent(inout)    :: grid
120!
121    type(Agrif_Grid_List), pointer :: sorted_child_list
122    type(Agrif_PGrid),     pointer :: child_p
123    type(Agrif_PGrid),     pointer :: g1p, g2p
124    type(Agrif_Proc_p),    pointer :: pp1, pp2
125    type(Agrif_Proc),      pointer :: proc
126    integer                        :: nb_seq_max, nb_seqs, cur_seq
127!
128    nb_seq_max = 0
129!
130    if ( grid % child_list % nitems == 0 ) return
131!
132! For each required proc...
133    pp1  => grid % required_proc_list % first
134    do while ( associated(pp1) )
135        proc => pp1 % proc
136        proc % nb_seqs = 0
137!   ..loop over all child grids...
138        child_p => grid % child_list % first
139        do while ( associated(child_p) )
140!       ..and look for 'proc' in the list of procs required by 'child'
141            pp2 => child_p % gr % required_proc_list % first
142            do while ( associated(pp2) )
143                if ( proc % pn == pp2 % proc % pn ) then
144!                   'proc' is required by this child grid, so we increment its number of sequences
145                    proc % nb_seqs = proc % nb_seqs + 1
146                    pp2 => NULL()
147                else
148                    pp2 => pp2 % next
149                endif
150            enddo
151            child_p => child_p % next
152        enddo
153        nb_seq_max = max(nb_seq_max, proc % nb_seqs)
154        pp1 => pp1 % next
155    enddo
156!
157! For each grid...
158    g1p => grid % child_list % first
159    do while ( associated(g1p) )
160!     compare it with the following ones
161        g2p => g1p % next
162        do while ( associated(g2p) )
163            if ( Agrif_seq_grids_are_connected( g1p % gr, g2p % gr ) ) then
164                call Agrif_gl_append( g1p % gr % neigh_list, g2p % gr )
165                call Agrif_gl_append( g2p % gr % neigh_list, g1p % gr )
166            endif
167            g2p => g2p % next
168        enddo
169        g1p => g1p % next
170    enddo
171!
172! Colorize graph nodes
173    nb_seqs = Agrif_seq_colorize_grid_list(grid % child_list)
174    sorted_child_list => Agrif_gl_merge_sort ( grid % child_list, compare_colors )
175!
176! Create sequence structure
177    cur_seq = 0
178    grid % child_seq => Agrif_seq_allocate_list(nb_seqs)
179    child_p => sorted_child_list % first
180    do while ( associated(child_p) )
181        if ( cur_seq /= child_p % gr % seq_num ) then
182            cur_seq = child_p % gr % seq_num
183        endif
184        call Agrif_seq_add_grid(grid % child_seq,cur_seq,child_p% gr)
185        child_p => child_p % next
186    enddo
187!
188    call Agrif_gl_delete(sorted_child_list)
189!---------------------------------------------------------------------------------------------------
190end subroutine Agrif_seq_create_proc_sequences
191!===================================================================================================
192!
193!===================================================================================================
194function Agrif_seq_grids_are_connected( g1, g2 ) result( connection )
195!
196!< Compare required_proc_list for g1 and g2. These are connected if they share a same proc.
197!---------------------------------------------------------------------------------------------------
198    type(Agrif_Grid), intent(in) :: g1, g2
199!
200    logical :: connection
201    type(Agrif_Proc_p), pointer :: pp1, pp2
202!
203    connection = .false.
204!
205    pp1 => g1 % required_proc_list % first
206!
207    do while( associated(pp1) .and. (.not. connection) )
208!
209        pp2 => g2 % required_proc_list % first
210        do while ( associated(pp2) .and. (.not. connection) )
211            if ( pp1 % proc % pn == pp2 % proc % pn ) then
212            ! if pp1 and pp2 are the same proc, it means that g1 and g2 are connected. We stop here.
213                connection = .true.
214            endif
215            pp2 => pp2 % next
216        enddo
217        pp1 => pp1 % next
218!
219    enddo
220!---------------------------------------------------------------------------------------------------
221end function Agrif_seq_grids_are_connected
222!===================================================================================================
223!
224!===================================================================================================
225function Agrif_seq_colorize_grid_list ( gridlist ) result ( ncolors )
226!
227!< 1.  Sort nodes in decreasing order of degree.
228!< 2.  Color the node with highest degree with color 1.
229!< 3.  Choose the node with the largest DSAT value. In case of conflict, choose the one with the
230!!       highest degree. Then the one corresponding to the largest grid.
231!< 4.  Color this node with the smallest possible color.
232!< 5.  If all nodes are colored, then stop. Otherwise, go to 3.
233!---------------------------------------------------------------------------------------------------
234    type(Agrif_Grid_List), intent(in)    :: gridlist
235!
236    type(Agrif_Grid_List),   pointer  :: X, Y
237    type(Agrif_PGrid),       pointer  :: gridp
238    type(Agrif_Grid_List),   pointer  :: tmp_gl
239    integer                           :: ncolors
240    integer, dimension(1:gridlist%nitems)   :: colors
241!
242! Be carefull...
243    nullify(Y)
244!
245! First initialize the color of each node
246    gridp => gridlist % first
247    do while ( associated(gridp) )
248        gridp % gr % seq_num = 0
249        gridp => gridp % next
250    enddo
251!
252! Then sort the grids by decreasing degree
253    X => Agrif_gl_merge_sort( gridlist, compare_grid_degrees )
254    gridp => X % first
255!
256! Colorize the first grid in the list
257    gridp % gr % seq_num = 1
258    gridp => gridp % next
259!
260! Then for each of its neighbors...
261    do while ( associated(gridp) )
262!
263        if ( gridp % gr % neigh_list % nitems == 0 ) then
264        ! this grid is alone... let.s attach it to an existing sequence
265            call Agrif_seq_attach_grid( X, gridp % gr )
266            gridp => gridp % next
267            cycle
268        endif
269!
270!     Compute dsat value of all non-colored grids
271        tmp_gl => Agrif_gl_build_from_gp(gridp)
272        call Agrif_seq_calc_dsat(tmp_gl)
273!
274!     Sort non-colored grids by decreasing dsat value, then by size
275        call Agrif_gl_delete(Y)
276        Y => Agrif_gl_merge_sort( tmp_gl, compare_dsat_values, compare_size_values )
277!
278!     Next coloration is for the first grid in this list  TODO : maybe we could find a better choice ..?
279        gridp => Y % first
280!
281!     Assign a color to the chosen grid
282        gridp % gr % seq_num = Agrif_seq_smallest_available_color_in_neighborhood( gridp % gr % neigh_list )
283!
284        gridp => gridp % next
285        call Agrif_gl_delete(tmp_gl)
286!
287    enddo
288!
289    call Agrif_gl_delete(X)
290    call Agrif_seq_colors_in_neighborhood( gridlist, colors )
291    ncolors = maxval(colors)
292!---------------------------------------------------------------------------------------------------
293end function Agrif_seq_colorize_grid_list
294!===================================================================================================
295!
296!===================================================================================================
297subroutine Agrif_seq_attach_grid ( gridlist, grid )
298!
299!< 'grid' is not connected to any neighbor. Therefore, we give an existing and well chosen color.
300!---------------------------------------------------------------------------------------------------
301    type(Agrif_Grid_List), intent(in)       :: gridlist
302    type(Agrif_Grid),      intent(inout)    :: grid
303!
304    integer, dimension(gridlist%nitems) :: colors
305    integer, dimension(:), allocatable  :: ngrids_by_color
306    integer :: i, color, ncolors
307!
308    call Agrif_seq_colors_in_neighborhood( gridlist, colors )
309    ncolors = maxval(colors)
310!
311    allocate(ngrids_by_color(ncolors))
312    ngrids_by_color = 0
313!
314    do i = 1,gridlist % nitems
315        if (colors(i) > 0)  ngrids_by_color(colors(i)) = ngrids_by_color(colors(i)) + 1
316    enddo
317!
318    color = ncolors
319    do i = 1,ncolors
320        if ( ngrids_by_color(i) < color ) color = i
321    enddo
322!
323    grid % seq_num = color
324    deallocate(ngrids_by_color)
325!---------------------------------------------------------------------------------------------------
326end subroutine Agrif_seq_attach_grid
327!===================================================================================================
328!
329!===================================================================================================
330subroutine Agrif_seq_colors_in_neighborhood ( neigh_list, colors )
331!---------------------------------------------------------------------------------------------------
332    type(Agrif_Grid_List), intent(in)   :: neigh_list
333    integer, dimension(:), intent(out)  :: colors
334!
335    integer                     :: i
336    type(Agrif_PGrid), pointer  :: gridp
337!
338    i = lbound(colors,1)
339    colors = 0
340    gridp => neigh_list % first
341!
342    do while ( associated(gridp) )
343!
344        if ( i > ubound(colors,1) ) then
345            print*,'Error in Agrif_seq_colors_in_neighborhood : "colors" array is too small'
346            stop
347        endif
348        colors(i) = gridp % gr % seq_num
349        gridp => gridp % next
350        i = i+1
351!
352    enddo
353!---------------------------------------------------------------------------------------------------
354end subroutine Agrif_seq_colors_in_neighborhood
355!===================================================================================================
356!
357!===================================================================================================
358function Agrif_seq_smallest_available_color_in_neighborhood ( neigh_list ) result ( smallest )
359!---------------------------------------------------------------------------------------------------
360    type(Agrif_Grid_List), intent(in)   :: neigh_list
361!
362    integer, dimension(:), allocatable :: color_is_met
363    integer     :: colors_tab(1:neigh_list%nitems)
364    integer     :: i, smallest, max_color
365!
366    call Agrif_seq_colors_in_neighborhood( neigh_list, colors_tab )
367    max_color = maxval(colors_tab)
368!
369    allocate(color_is_met(1:max_color))
370    color_is_met = 0
371!
372    do i = 1,neigh_list % nitems
373        if ( colors_tab(i) /= 0 ) then
374            color_is_met(colors_tab(i)) = 1
375        endif
376    enddo
377!
378    smallest = max_color+1
379    do i = 1,max_color
380        if ( color_is_met(i) == 0 ) then
381            smallest = i
382            exit
383        endif
384    enddo
385!
386    deallocate(color_is_met)
387!---------------------------------------------------------------------------------------------------
388end function Agrif_seq_smallest_available_color_in_neighborhood
389!===================================================================================================
390!
391!===================================================================================================
392subroutine Agrif_seq_calc_dsat ( gridlist )
393!< For each node 'v' :
394!<   if none of its neighbors is colored then
395!<       DSAT(v) = degree(v)  #  degree(v) := number of neighbors
396!<   else
397!<       DSAT(v) = number of different colors used in the first neighborhood of v.
398!---------------------------------------------------------------------------------------------------
399    type(Agrif_Grid_List), intent(in)   :: gridlist
400!
401    type(Agrif_PGrid), pointer          :: gridp
402    type(Agrif_Grid),  pointer          :: grid
403    integer, dimension(:), allocatable  :: colors, color_is_met
404    integer                             :: i, ncolors
405!
406    gridp => gridlist % first
407!
408    do while ( associated(gridp) )
409!
410        grid => gridp % gr
411!
412        allocate(colors(grid % neigh_list % nitems))
413        call Agrif_seq_colors_in_neighborhood( grid % neigh_list, colors )
414
415        allocate(color_is_met(1:maxval(colors)))
416        color_is_met = 0
417!
418        do i = 1,grid % neigh_list % nitems
419            if ( colors(i) /= 0 ) then
420                color_is_met(colors(i)) = 1
421            endif
422        enddo
423        ncolors = sum(color_is_met)
424!
425        if ( ncolors == 0 ) then
426            grid % dsat = grid % neigh_list % nitems
427        else
428            grid % dsat = ncolors
429        endif
430        deallocate(colors, color_is_met)
431        gridp => gridp % next
432    enddo
433!---------------------------------------------------------------------------------------------------
434end subroutine Agrif_seq_calc_dsat
435!===================================================================================================
436!
437!===================================================================================================
438subroutine Agrif_seq_allocate_procs_to_childs ( coarse_grid )
439!---------------------------------------------------------------------------------------------------
440    type(Agrif_Grid), intent(inout) :: coarse_grid
441!
442    integer                         :: is, ip, ig, ngrids
443    type(Agrif_Grid_List), pointer  :: gridlist
444    type(Agrif_PGrid),     pointer  :: gp
445    type(Agrif_Grid),      pointer  :: grid
446    type(Agrif_Proc_List), pointer  :: proclist
447    type(Agrif_Proc),      pointer  :: proc
448    type(Agrif_Proc_p),    pointer  :: pp
449    type(Agrif_Proc), dimension(:), allocatable, target :: procarray
450    type(Agrif_Grid), dimension(:), allocatable         :: gridarray
451    type(Agrif_Sequence_List), pointer :: seqlist
452    real,dimension(:),allocatable :: grid_costs
453    integer,dimension(:), allocatable :: nbprocs_per_grid
454    integer :: i1, i2, j1, j2
455    real :: max_cost
456    integer :: max_index
457!
458    seqlist => coarse_grid % child_seq
459    if ( .not. associated(seqlist) ) return
460!
461! Initialize proc allocation
462    pp => coarse_grid % proc_def_list % first
463    do while ( associated(pp) )
464        pp % proc % grid_id = 0
465        pp => pp % next
466    enddo
467!
468! For each sequence...
469    do is = 1,seqlist % nb_seqs
470!
471        proclist => seqlist % sequences(is) % proclist
472        gridlist => seqlist % sequences(is) % gridlist
473!
474!     Copy coarse grid proc list and convert it to an array
475        call Agrif_pl_deep_copy( coarse_grid % proc_def_list, proclist )
476        call Agrif_pl_to_array ( proclist, procarray )
477!
478!     Allocate a temporary array with concerned grid numbers
479        ngrids = gridlist % nitems
480        allocate(gridarray(1:ngrids))
481        allocate(grid_costs(1:ngrids))
482        allocate(nbprocs_per_grid(1:ngrids))
483
484      nbprocs_per_grid = 0
485!
486!     Allocate required procs to each grid
487        gp => gridlist % first
488        ig = 0
489        do while ( associated(gp) )
490            grid => gp % gr
491            ig = ig+1 ; gridarray(ig) = grid
492            pp => grid % required_proc_list % first
493            do while ( associated(pp) )
494                procarray( pp % proc % pn+1 ) % grid_id = grid % fixedrank
495                nbprocs_per_grid(ig) = nbprocs_per_grid(ig) + 1
496                pp => pp % next
497            enddo
498            gp => gp % next
499        enddo
500!
501!     Add unused procs to the grids
502! TODO FIXME: This is just a dummy allocation. You should take into account grid size and more
503!             information here...
504
505! Estimate current costs
506
507        do ig = 1, ngrids
508          i1 = gridarray(ig)%ix(1)
509          i2 = gridarray(ig)%ix(1)+gridarray(ig)%nb(1)/gridarray(ig)%spaceref(1)-1
510          j1 = gridarray(ig)%ix(2)
511          j2 = gridarray(ig)%ix(2)+gridarray(ig)%nb(2)/gridarray(ig)%spaceref(2)-1
512          Call Agrif_estimate_parallel_cost(i1,i2,j1,j2,nbprocs_per_grid(ig),grid_costs(ig))
513          grid_costs(ig) = grid_costs(ig) * gridarray(ig)%timeref(1)
514        enddo
515
516        ig = 1
517        do ip = 1,proclist%nitems
518            if ( procarray( ip ) % grid_id == 0 ) then
519!             this proc is unused
520              max_cost = 0.
521              max_index = 1
522              do ig = 1,ngrids
523                 if (grid_costs(ig) > max_cost) then
524                    max_cost = grid_costs(ig)
525                    max_index = ig
526                 endif
527              enddo
528
529              ig = max_index
530              procarray( ip ) % grid_id = gridarray(ig) % fixedrank
531
532              nbprocs_per_grid(ig) =  nbprocs_per_grid(ig) + 1
533              i1 = gridarray(ig)%ix(1)
534              i2 = gridarray(ig)%ix(1)+gridarray(ig)%nb(1)/gridarray(ig)%spaceref(1)-1
535              j1 = gridarray(ig)%ix(2)
536              j2 = gridarray(ig)%ix(2)+gridarray(ig)%nb(2)/gridarray(ig)%spaceref(2)-1
537              Call Agrif_estimate_parallel_cost(i1,i2,j1,j2,nbprocs_per_grid(ig),grid_costs(ig))
538              grid_costs(ig) = grid_costs(ig) * gridarray(ig)%timeref(1)
539
540            endif
541        enddo
542!
543!     Allocate proc nums to each grid
544        gp => gridlist % first
545        do while ( associated(gp) )
546            do ip = 1,proclist%nitems
547                if ( procarray( ip ) % grid_id == gp % gr % fixedrank ) then
548                    allocate(proc)
549                    proc = procarray( ip )
550                    call Agrif_pl_append(gp % gr % proc_def_in_parent_list, proc)
551                endif
552            enddo
553            gp => gp % next
554        enddo
555!
556!     Clean up
557        deallocate(procarray, gridarray, grid_costs, nbprocs_per_grid)
558!
559    enddo
560!---------------------------------------------------------------------------------------------------
561end subroutine Agrif_seq_allocate_procs_to_childs
562!===================================================================================================
563!
564!===================================================================================================
565subroutine Agrif_seq_create_communicators ( grid )
566!---------------------------------------------------------------------------------------------------
567    type(Agrif_Grid), intent(inout) :: grid
568!
569    include 'mpif.h'
570    type(Agrif_Sequence_List), pointer  :: seqlist  ! List of child sequences
571    type(Agrif_PGrid), pointer          :: gridp
572    type(Agrif_Proc),  pointer          :: proc
573    integer     :: i, ierr
574    integer     :: current_comm, comm_seq, color_seq
575!
576    seqlist => grid % child_seq
577    if ( .not. associated(seqlist) ) return
578!
579    current_comm = grid % communicator
580    color_seq = MPI_COMM_NULL
581!
582! For each sequence, split the current communicator into as many groups as needed.
583    do i = 1,seqlist % nb_seqs
584!
585!     Loop over each proclist to give a color to the current process
586        gridp => seqlist % sequences(i) % gridlist % first
587        grid_loop : do while ( associated(gridp) )
588            proc => Agrif_pl_search_proc( gridp % gr % proc_def_in_parent_list, Agrif_Procrank )
589            if ( associated(proc) ) then
590                if ( gridp % gr % fixedrank /= proc % grid_id ) then
591                    write(*,'("### Error Agrif_seq_create_communicators : ")')
592                    write(*,'("  inconsitancy on proc ",i2,":")') Agrif_Procrank
593                    write(*,'("gr % fixedrank = ",i0,", where proc % grid_id = ",i0)') &
594                            gridp % gr % fixedrank, proc % grid_id
595                    stop
596                endif
597                color_seq = gridp % gr % fixedrank
598                exit grid_loop
599            endif
600            gridp => gridp % next
601        enddo grid_loop
602!
603        call MPI_COMM_SPLIT(current_comm, color_seq, Agrif_ProcRank, comm_seq, ierr)
604        gridp % gr % communicator = comm_seq
605!
606    enddo
607!---------------------------------------------------------------------------------------------------
608end subroutine Agrif_seq_create_communicators
609!===================================================================================================
610!
611!===================================================================================================
612function Agrif_seq_select_child ( g, is ) result ( gridp )
613!---------------------------------------------------------------------------------------------------
614    type(Agrif_Grid), pointer, intent(in)   :: g
615    integer,                   intent(in)   :: is
616!
617    type(Agrif_PGrid), pointer  :: gridp
618    type(Agrif_Proc),  pointer  :: proc
619!
620    call Agrif_Instance( g )
621    gridp => g % child_seq % sequences(is) % gridlist % first
622!
623    do while ( associated(gridp) )
624        proc => Agrif_pl_search_proc( gridp % gr % proc_def_in_parent_list, Agrif_Procrank )
625        if ( associated(proc) ) then
626            return
627        endif
628        gridp => gridp % next
629    enddo
630    write(*,'("### Error Agrif_seq_select_child : no grid found in sequence ",i0," (mother G",i0,") for P",i0)')&
631        is, g%fixedrank, Agrif_Procrank
632    stop
633!---------------------------------------------------------------------------------------------------
634end function Agrif_seq_select_child
635!===================================================================================================
636#else
637    subroutine dummy_Agrif_seq ()
638    end subroutine dummy_Agrif_seq
639#endif
640!
641end module Agrif_seq
Note: See TracBrowser for help on using the repository browser.