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