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 | ! |
---|
24 | !> Module Agrif_Save |
---|
25 | !! |
---|
26 | !! Module for the initialization by copy of the grids created by clustering. |
---|
27 | ! |
---|
28 | module Agrif_Save |
---|
29 | ! |
---|
30 | use Agrif_Types |
---|
31 | use Agrif_Link |
---|
32 | use Agrif_Arrays |
---|
33 | use Agrif_Variables |
---|
34 | ! |
---|
35 | implicit none |
---|
36 | ! |
---|
37 | contains |
---|
38 | ! |
---|
39 | !=================================================================================================== |
---|
40 | ! subroutine Agrif_deallocate_Arrays |
---|
41 | !--------------------------------------------------------------------------------------------------- |
---|
42 | subroutine Agrif_deallocate_Arrays ( var ) |
---|
43 | !--------------------------------------------------------------------------------------------------- |
---|
44 | type(Agrif_Variable), pointer :: var |
---|
45 | ! |
---|
46 | if (allocated(var%array1)) deallocate(var%array1) |
---|
47 | if (allocated(var%array2)) deallocate(var%array2) |
---|
48 | if (allocated(var%array3)) deallocate(var%array3) |
---|
49 | if (allocated(var%array4)) deallocate(var%array4) |
---|
50 | if (allocated(var%array5)) deallocate(var%array5) |
---|
51 | if (allocated(var%array6)) deallocate(var%array6) |
---|
52 | ! |
---|
53 | if (allocated(var%darray1)) deallocate(var%darray1) |
---|
54 | if (allocated(var%darray2)) deallocate(var%darray2) |
---|
55 | if (allocated(var%darray3)) deallocate(var%darray3) |
---|
56 | if (allocated(var%darray4)) deallocate(var%darray4) |
---|
57 | if (allocated(var%darray5)) deallocate(var%darray5) |
---|
58 | if (allocated(var%darray6)) deallocate(var%darray6) |
---|
59 | ! |
---|
60 | if (allocated(var%sarray1)) deallocate(var%sarray1) |
---|
61 | if (allocated(var%sarray2)) deallocate(var%sarray2) |
---|
62 | if (allocated(var%sarray3)) deallocate(var%sarray3) |
---|
63 | if (allocated(var%sarray4)) deallocate(var%sarray4) |
---|
64 | if (allocated(var%sarray5)) deallocate(var%sarray5) |
---|
65 | if (allocated(var%sarray6)) deallocate(var%sarray6) |
---|
66 | ! |
---|
67 | ! |
---|
68 | if (associated(var%oldvalues2D)) deallocate(var%oldvalues2D) |
---|
69 | if (allocated(var%posvar)) deallocate(var%posvar) |
---|
70 | if (allocated(var%interptab)) deallocate(var%interptab) |
---|
71 | if (allocated(var%coords)) deallocate(var%coords) |
---|
72 | !--------------------------------------------------------------------------------------------------- |
---|
73 | end subroutine Agrif_deallocate_Arrays |
---|
74 | !--------------------------------------------------------------------------------------------------- |
---|
75 | subroutine Agrif_deallocate_Arrays_c ( var_c ) |
---|
76 | !--------------------------------------------------------------------------------------------------- |
---|
77 | type(Agrif_Variable_c), pointer :: var_c |
---|
78 | ! |
---|
79 | if (allocated(var_c%carray1)) deallocate(var_c%carray1) |
---|
80 | if (allocated(var_C%carray2)) deallocate(var_c%carray2) |
---|
81 | if (allocated(var_C%carrayu)) deallocate(var_c%carrayu) |
---|
82 | ! |
---|
83 | !--------------------------------------------------------------------------------------------------- |
---|
84 | end subroutine Agrif_deallocate_Arrays_c |
---|
85 | !=================================================================================================== |
---|
86 | !--------------------------------------------------------------------------------------------------- |
---|
87 | subroutine Agrif_deallocate_Arrays_l ( var_l ) |
---|
88 | !--------------------------------------------------------------------------------------------------- |
---|
89 | type(Agrif_Variable_l), pointer :: var_l |
---|
90 | ! |
---|
91 | ! |
---|
92 | if (allocated(var_l%larray1)) deallocate(var_l%larray1) |
---|
93 | if (allocated(var_l%larray2)) deallocate(var_l%larray2) |
---|
94 | if (allocated(var_l%larray3)) deallocate(var_l%larray3) |
---|
95 | if (allocated(var_l%larray4)) deallocate(var_l%larray4) |
---|
96 | if (allocated(var_l%larray5)) deallocate(var_l%larray5) |
---|
97 | if (allocated(var_l%larray6)) deallocate(var_l%larray6) |
---|
98 | ! |
---|
99 | !--------------------------------------------------------------------------------------------------- |
---|
100 | end subroutine Agrif_deallocate_Arrays_l |
---|
101 | !--------------------------------------------------------------------------------------------------- |
---|
102 | subroutine Agrif_deallocate_Arrays_i ( var_i ) |
---|
103 | !--------------------------------------------------------------------------------------------------- |
---|
104 | type(Agrif_Variable_i), pointer :: var_i |
---|
105 | ! |
---|
106 | ! |
---|
107 | if (allocated(var_i%iarray1)) deallocate(var_i%iarray1) |
---|
108 | if (allocated(var_i%iarray2)) deallocate(var_i%iarray2) |
---|
109 | if (allocated(var_i%iarray3)) deallocate(var_i%iarray3) |
---|
110 | if (allocated(var_i%iarray4)) deallocate(var_i%iarray4) |
---|
111 | if (allocated(var_i%iarray5)) deallocate(var_i%iarray5) |
---|
112 | if (allocated(var_i%iarray6)) deallocate(var_i%iarray6) |
---|
113 | ! |
---|
114 | !--------------------------------------------------------------------------------------------------- |
---|
115 | end subroutine Agrif_deallocate_Arrays_i |
---|
116 | ! |
---|
117 | !=================================================================================================== |
---|
118 | ! subroutine Agrif_Free_data_before |
---|
119 | !--------------------------------------------------------------------------------------------------- |
---|
120 | subroutine Agrif_Free_data_before ( Agrif_Gr ) |
---|
121 | !--------------------------------------------------------------------------------------------------- |
---|
122 | type(Agrif_Grid), pointer :: Agrif_Gr ! Pointer on the current grid |
---|
123 | ! |
---|
124 | integer :: i |
---|
125 | type(Agrif_Variables_List), pointer :: parcours |
---|
126 | type(Agrif_Variable), pointer :: var |
---|
127 | type(Agrif_Variable_c), pointer :: var_c |
---|
128 | type(Agrif_Variable_r), pointer :: var_r |
---|
129 | type(Agrif_Variable_l), pointer :: var_l |
---|
130 | type(Agrif_Variable_i), pointer :: var_i |
---|
131 | ! |
---|
132 | do i = 1,Agrif_NbVariables(0) |
---|
133 | ! |
---|
134 | var => Agrif_Gr % tabvars(i) |
---|
135 | ! |
---|
136 | if ( .NOT. Agrif_Mygrid % tabvars(i) % restore ) then |
---|
137 | call Agrif_deallocate_Arrays(var) |
---|
138 | endif |
---|
139 | ! |
---|
140 | if (associated(var%list_interp)) then |
---|
141 | call Agrif_Free_list_interp(var%list_interp) |
---|
142 | endif |
---|
143 | ! |
---|
144 | enddo |
---|
145 | do i=1,Agrif_NbVariables(1) |
---|
146 | var_c => Agrif_Gr % tabvars_c(i) |
---|
147 | call Agrif_deallocate_Arrays_c(var_c) |
---|
148 | enddo |
---|
149 | do i=1,Agrif_NbVariables(3) |
---|
150 | var_l => Agrif_Gr % tabvars_l(i) |
---|
151 | call Agrif_deallocate_Arrays_l(var_l) |
---|
152 | enddo |
---|
153 | do i=1,Agrif_NbVariables(4) |
---|
154 | var_i => Agrif_Gr % tabvars_i(i) |
---|
155 | call Agrif_deallocate_Arrays_i(var_i) |
---|
156 | enddo |
---|
157 | |
---|
158 | parcours => Agrif_Gr % variables |
---|
159 | |
---|
160 | do i = 1,Agrif_Gr%NbVariables |
---|
161 | ! |
---|
162 | if ( .NOT. parcours%var%root_var % restore ) then |
---|
163 | call Agrif_deallocate_Arrays(parcours%var) |
---|
164 | endif |
---|
165 | ! |
---|
166 | if (associated(parcours%var%list_interp)) then |
---|
167 | call Agrif_Free_list_interp(parcours%var%list_interp) |
---|
168 | endif |
---|
169 | ! |
---|
170 | if ( .NOT. parcours%var%root_var % restore ) then |
---|
171 | deallocate(parcours%var) |
---|
172 | endif |
---|
173 | ! |
---|
174 | parcours => parcours % next |
---|
175 | ! |
---|
176 | enddo |
---|
177 | ! |
---|
178 | if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) then |
---|
179 | if ( Agrif_Probdim == 1 ) deallocate(Agrif_Gr%tabpoint1D) |
---|
180 | if ( Agrif_Probdim == 2 ) deallocate(Agrif_Gr%tabpoint2D) |
---|
181 | if ( Agrif_Probdim == 3 ) deallocate(Agrif_Gr%tabpoint3D) |
---|
182 | endif |
---|
183 | !--------------------------------------------------------------------------------------------------- |
---|
184 | end subroutine Agrif_Free_data_before |
---|
185 | !=================================================================================================== |
---|
186 | ! |
---|
187 | !=================================================================================================== |
---|
188 | ! subroutine Agrif_Free_list_interp |
---|
189 | !--------------------------------------------------------------------------------------------------- |
---|
190 | recursive subroutine Agrif_Free_list_interp ( list_interp ) |
---|
191 | !--------------------------------------------------------------------------------------------------- |
---|
192 | type(Agrif_List_Interp_Loc), pointer :: list_interp |
---|
193 | ! |
---|
194 | if (associated(list_interp%suiv)) call Agrif_Free_list_interp(list_interp%suiv) |
---|
195 | |
---|
196 | #if defined AGRIF_MPI |
---|
197 | deallocate(list_interp%interp_loc%tab4t) |
---|
198 | deallocate(list_interp%interp_loc%memberinall) |
---|
199 | deallocate(list_interp%interp_loc%sendtoproc1) |
---|
200 | deallocate(list_interp%interp_loc%recvfromproc1) |
---|
201 | #endif |
---|
202 | deallocate(list_interp%interp_loc) |
---|
203 | deallocate(list_interp) |
---|
204 | nullify(list_interp) |
---|
205 | !--------------------------------------------------------------------------------------------------- |
---|
206 | end subroutine Agrif_Free_list_interp |
---|
207 | !=================================================================================================== |
---|
208 | ! |
---|
209 | !=================================================================================================== |
---|
210 | ! subroutine Agrif_Free_data_after |
---|
211 | !--------------------------------------------------------------------------------------------------- |
---|
212 | subroutine Agrif_Free_data_after ( Agrif_Gr ) |
---|
213 | !--------------------------------------------------------------------------------------------------- |
---|
214 | type(Agrif_Grid), pointer :: Agrif_Gr !< Pointer on the current grid |
---|
215 | ! |
---|
216 | integer :: i |
---|
217 | type(Agrif_Variables_List), pointer :: parcours, rootparcours |
---|
218 | type(Agrif_Variable), pointer :: var |
---|
219 | ! |
---|
220 | do i = 1,Agrif_NbVariables(0) |
---|
221 | if ( Agrif_Mygrid % tabvars(i) % restore ) then |
---|
222 | var => Agrif_Gr%tabvars(i) |
---|
223 | call Agrif_deallocate_Arrays(var) |
---|
224 | endif |
---|
225 | enddo |
---|
226 | ! |
---|
227 | parcours => Agrif_Gr%variables |
---|
228 | rootparcours => Agrif_Mygrid%variables |
---|
229 | ! |
---|
230 | do i = 1,Agrif_Gr%NbVariables |
---|
231 | if (rootparcours%var % restore ) then |
---|
232 | call Agrif_deallocate_Arrays(parcours%var) |
---|
233 | deallocate(parcours%var) |
---|
234 | endif |
---|
235 | parcours => parcours % next |
---|
236 | rootparcours => rootparcours % next |
---|
237 | enddo |
---|
238 | !--------------------------------------------------------------------------------------------------- |
---|
239 | end subroutine Agrif_Free_data_after |
---|
240 | !=================================================================================================== |
---|
241 | ! |
---|
242 | !=================================================================================================== |
---|
243 | ! subroutine Agrif_CopyFromOld_All |
---|
244 | ! |
---|
245 | !> Called in Agrif_Clustering#Agrif_Init_Hierarchy. |
---|
246 | !--------------------------------------------------------------------------------------------------- |
---|
247 | recursive subroutine Agrif_CopyFromOld_All ( g, oldchildlist ) |
---|
248 | !--------------------------------------------------------------------------------------------------- |
---|
249 | type(Agrif_Grid), pointer, intent(inout) :: g !< Pointer on the current grid |
---|
250 | type(Agrif_Grid_List), intent(in) :: oldchildlist |
---|
251 | ! |
---|
252 | type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure |
---|
253 | real :: g_eps, eps, oldgrid_eps |
---|
254 | integer :: out |
---|
255 | integer :: iii |
---|
256 | ! |
---|
257 | out = 0 |
---|
258 | ! |
---|
259 | parcours => oldchildlist % first |
---|
260 | ! |
---|
261 | do while (associated(parcours)) |
---|
262 | ! |
---|
263 | if ((.NOT. g % fixed) .AND. (parcours % gr %oldgrid)) then |
---|
264 | ! |
---|
265 | g_eps = huge(1.) |
---|
266 | oldgrid_eps = huge(1.) |
---|
267 | do iii = 1,Agrif_Probdim |
---|
268 | g_eps = min(g_eps,g % Agrif_dx(iii)) |
---|
269 | oldgrid_eps = min(oldgrid_eps, parcours % gr % Agrif_dx(iii)) |
---|
270 | enddo |
---|
271 | ! |
---|
272 | eps = min(g_eps,oldgrid_eps)/100. |
---|
273 | ! |
---|
274 | do iii = 1,Agrif_Probdim |
---|
275 | if (g % Agrif_dx(iii) < (parcours % gr % Agrif_dx(iii) - eps)) then |
---|
276 | parcours => parcours % next |
---|
277 | out = 1 |
---|
278 | exit |
---|
279 | endif |
---|
280 | enddo |
---|
281 | ! |
---|
282 | if ( out == 1 ) cycle |
---|
283 | ! |
---|
284 | call Agrif_CopyFromOld(g,parcours%gr) |
---|
285 | ! |
---|
286 | endif |
---|
287 | ! |
---|
288 | call Agrif_CopyFromOld_All(g,parcours % gr % child_list) |
---|
289 | ! |
---|
290 | parcours => parcours % next |
---|
291 | ! |
---|
292 | enddo |
---|
293 | !--------------------------------------------------------------------------------------------------- |
---|
294 | end subroutine Agrif_CopyFromOld_All |
---|
295 | !=================================================================================================== |
---|
296 | ! |
---|
297 | !=================================================================================================== |
---|
298 | ! subroutine Agrif_CopyFromOld |
---|
299 | ! |
---|
300 | !> Call to the Agrif_Copy procedure. |
---|
301 | !--------------------------------------------------------------------------------------------------- |
---|
302 | subroutine Agrif_CopyFromOld ( new_gr, old_gr ) |
---|
303 | !--------------------------------------------------------------------------------------------------- |
---|
304 | type(Agrif_Grid), pointer, intent(inout) :: new_gr !< Pointer on the current grid |
---|
305 | type(Agrif_Grid), pointer, intent(inout) :: old_gr !< Pointer on an old grid |
---|
306 | ! |
---|
307 | type(Agrif_Variable), pointer :: new_var |
---|
308 | type(Agrif_Variable), pointer :: old_var |
---|
309 | integer :: i |
---|
310 | ! |
---|
311 | do i = 1,Agrif_NbVariables(0) |
---|
312 | if ( Agrif_Mygrid % tabvars(i) % restore ) then |
---|
313 | old_var => old_gr % tabvars(i) |
---|
314 | new_var => new_gr % tabvars(i) |
---|
315 | call Agrif_Copy(new_gr, old_gr, new_var, old_var ) |
---|
316 | endif |
---|
317 | enddo |
---|
318 | !--------------------------------------------------------------------------------------------------- |
---|
319 | end subroutine Agrif_CopyFromOld |
---|
320 | !=================================================================================================== |
---|
321 | ! |
---|
322 | !=================================================================================================== |
---|
323 | ! subroutine Agrif_CopyFromOld_AllOneVar |
---|
324 | ! |
---|
325 | !> Called in Agrif_Clustering#Agrif_Init_Hierarchy. |
---|
326 | !--------------------------------------------------------------------------------------------------- |
---|
327 | recursive subroutine Agrif_CopyFromOld_AllOneVar ( g, oldchildlist, indic ) |
---|
328 | !--------------------------------------------------------------------------------------------------- |
---|
329 | type(Agrif_Grid), pointer, intent(inout) :: g !< Pointer on the current grid |
---|
330 | type(Agrif_Grid_List), intent(in) :: oldchildlist |
---|
331 | integer, intent(in) :: indic |
---|
332 | ! |
---|
333 | type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure |
---|
334 | real :: g_eps,eps,oldgrid_eps |
---|
335 | integer :: out |
---|
336 | integer :: iii |
---|
337 | ! |
---|
338 | out = 0 |
---|
339 | ! |
---|
340 | parcours => oldchildlist % first |
---|
341 | ! |
---|
342 | do while (associated(parcours)) |
---|
343 | ! |
---|
344 | if ((.NOT. g % fixed) .AND. (parcours % gr %oldgrid)) then |
---|
345 | ! |
---|
346 | g_eps = huge(1.) |
---|
347 | oldgrid_eps = huge(1.) |
---|
348 | do iii = 1,Agrif_Probdim |
---|
349 | g_eps = min(g_eps,g % Agrif_dx(iii)) |
---|
350 | oldgrid_eps = min(oldgrid_eps,parcours % gr % Agrif_dx(iii)) |
---|
351 | enddo |
---|
352 | ! |
---|
353 | eps = min(g_eps,oldgrid_eps)/100. |
---|
354 | ! |
---|
355 | do iii = 1,Agrif_Probdim |
---|
356 | if (g % Agrif_dx(iii) < (parcours % gr % Agrif_dx(iii) - eps)) then |
---|
357 | parcours => parcours % next |
---|
358 | out = 1 |
---|
359 | exit |
---|
360 | endif |
---|
361 | enddo |
---|
362 | |
---|
363 | if ( out == 1 ) cycle |
---|
364 | ! |
---|
365 | call Agrif_CopyFromOldOneVar(g,parcours%gr,indic) |
---|
366 | ! |
---|
367 | endif |
---|
368 | ! |
---|
369 | call Agrif_CopyFromOld_AllOneVar(g,parcours%gr % child_list,indic) |
---|
370 | ! |
---|
371 | parcours => parcours % next |
---|
372 | ! |
---|
373 | enddo |
---|
374 | !--------------------------------------------------------------------------------------------------- |
---|
375 | end subroutine Agrif_CopyFromOld_AllOneVar |
---|
376 | !=================================================================================================== |
---|
377 | ! |
---|
378 | !=================================================================================================== |
---|
379 | ! subroutine Agrif_CopyFromOldOneVar |
---|
380 | ! |
---|
381 | !> Call to Agrif_Copy |
---|
382 | !--------------------------------------------------------------------------------------------------- |
---|
383 | subroutine Agrif_CopyFromOldOneVar ( new_gr, old_gr, indic ) |
---|
384 | !--------------------------------------------------------------------------------------------------- |
---|
385 | type(Agrif_Grid), pointer, intent(inout) :: new_gr !< Pointer on the current grid |
---|
386 | type(Agrif_Grid), pointer, intent(in) :: old_gr !< Pointer on an old grid |
---|
387 | integer, intent(in) :: indic |
---|
388 | ! |
---|
389 | type(Agrif_Variable), pointer :: new_var, old_var |
---|
390 | ! |
---|
391 | new_var => Agrif_Search_Variable(new_gr, -indic) |
---|
392 | old_var => Agrif_Search_Variable(old_gr, -indic) |
---|
393 | ! |
---|
394 | call Agrif_array_allocate(new_var,new_var%lb,new_var%ub,new_var%nbdim) |
---|
395 | call Agrif_Copy(new_gr, old_gr, new_var,old_var) |
---|
396 | !--------------------------------------------------------------------------------------------------- |
---|
397 | end subroutine Agrif_CopyFromOldOneVar |
---|
398 | !=================================================================================================== |
---|
399 | ! |
---|
400 | !=================================================================================================== |
---|
401 | ! subroutine Agrif_Copy |
---|
402 | !--------------------------------------------------------------------------------------------------- |
---|
403 | subroutine Agrif_Copy ( new_gr, old_gr, new_var, old_var ) |
---|
404 | !--------------------------------------------------------------------------------------------------- |
---|
405 | type(Agrif_Grid), pointer, intent(in) :: new_gr !< Pointer on the current grid |
---|
406 | type(Agrif_Grid), pointer, intent(in) :: old_gr !< Pointer on an old grid |
---|
407 | type(Agrif_Variable), pointer, intent(inout) :: new_var |
---|
408 | type(Agrif_Variable), pointer, intent(in) :: old_var |
---|
409 | ! |
---|
410 | type(Agrif_Variable), pointer :: root ! Pointer on the variable of the root grid |
---|
411 | integer :: nbdim ! Number of dimensions of the current grid |
---|
412 | integer, dimension(6) :: pttabnew ! Indexes of the first point in the domain |
---|
413 | integer, dimension(6) :: petabnew ! Indexes of the first point in the domain |
---|
414 | integer, dimension(6) :: pttabold ! Indexes of the first point in the domain |
---|
415 | integer, dimension(6) :: petabold ! Indexes of the first point in the domain |
---|
416 | integer, dimension(6) :: nbtabold ! Number of cells in each direction |
---|
417 | integer, dimension(6) :: nbtabnew ! Number of cells in each direction |
---|
418 | real, dimension(6) :: snew,sold |
---|
419 | real, dimension(6) :: dsnew,dsold |
---|
420 | real :: eps |
---|
421 | integer :: n |
---|
422 | ! |
---|
423 | root => new_var % root_var |
---|
424 | nbdim = root % nbdim |
---|
425 | ! |
---|
426 | do n = 1,nbdim |
---|
427 | ! |
---|
428 | select case(root % interptab(n)) |
---|
429 | ! |
---|
430 | case('x') |
---|
431 | ! |
---|
432 | pttabnew(n) = root % point(n) |
---|
433 | pttabold(n) = root % point(n) |
---|
434 | snew(n) = new_gr % Agrif_x(1) |
---|
435 | sold(n) = old_gr % Agrif_x(1) |
---|
436 | dsnew(n) = new_gr % Agrif_dx(1) |
---|
437 | dsold(n) = old_gr % Agrif_dx(1) |
---|
438 | ! |
---|
439 | if (root % posvar(n) == 1) then |
---|
440 | petabnew(n) = pttabnew(n) + new_gr % nb(1) |
---|
441 | petabold(n) = pttabold(n) + old_gr % nb(1) |
---|
442 | else |
---|
443 | petabnew(n) = pttabnew(n) + new_gr % nb(1) - 1 |
---|
444 | petabold(n) = pttabold(n) + old_gr % nb(1) - 1 |
---|
445 | snew(n) = snew(n) + dsnew(n)/2. |
---|
446 | sold(n) = sold(n) + dsold(n)/2. |
---|
447 | endif |
---|
448 | ! |
---|
449 | case('y') |
---|
450 | ! |
---|
451 | pttabnew(n) = root % point(n) |
---|
452 | pttabold(n) = root % point(n) |
---|
453 | snew(n) = new_gr % Agrif_x(2) |
---|
454 | sold(n) = old_gr % Agrif_x(2) |
---|
455 | dsnew(n) = new_gr % Agrif_dx(2) |
---|
456 | dsold(n) = old_gr % Agrif_dx(2) |
---|
457 | ! |
---|
458 | if (root % posvar(n) == 1) then |
---|
459 | petabnew(n) = pttabnew(n) + new_gr % nb(2) |
---|
460 | petabold(n) = pttabold(n) + old_gr % nb(2) |
---|
461 | else |
---|
462 | petabnew(n) = pttabnew(n) + new_gr % nb(2) - 1 |
---|
463 | petabold(n) = pttabold(n) + old_gr % nb(2) - 1 |
---|
464 | snew(n) = snew(n) + dsnew(n)/2. |
---|
465 | sold(n) = sold(n) + dsold(n)/2. |
---|
466 | endif |
---|
467 | ! |
---|
468 | case('z') |
---|
469 | ! |
---|
470 | pttabnew(n) = root % point(n) |
---|
471 | pttabold(n) = root % point(n) |
---|
472 | snew(n) = new_gr % Agrif_x(3) |
---|
473 | sold(n) = old_gr % Agrif_x(3) |
---|
474 | dsnew(n) = new_gr % Agrif_dx(3) |
---|
475 | dsold(n) = old_gr % Agrif_dx(3) |
---|
476 | ! |
---|
477 | if (root % posvar(n) == 1) then |
---|
478 | petabnew(n) = pttabnew(n) + new_gr % nb(3) |
---|
479 | petabold(n) = pttabold(n) + old_gr % nb(3) |
---|
480 | else |
---|
481 | petabnew(n) = pttabnew(n) + new_gr % nb(3) - 1 |
---|
482 | petabold(n) = pttabold(n) + old_gr % nb(3) - 1 |
---|
483 | snew(n) = snew(n) + dsnew(n)/2. |
---|
484 | sold(n) = sold(n) + dsold(n)/2. |
---|
485 | endif |
---|
486 | ! |
---|
487 | case('N') |
---|
488 | ! |
---|
489 | call Agrif_get_var_bounds(new_var,pttabnew(n),petabnew(n),n) |
---|
490 | ! |
---|
491 | pttabold(n) = pttabnew(n) |
---|
492 | petabold(n) = petabnew(n) |
---|
493 | snew(n) = 0. |
---|
494 | sold(n) = 0. |
---|
495 | dsnew(n) = 1. |
---|
496 | dsold(n) = 1. |
---|
497 | ! |
---|
498 | end select |
---|
499 | ! |
---|
500 | enddo |
---|
501 | ! |
---|
502 | do n = 1,nbdim |
---|
503 | nbtabnew(n) = petabnew(n) - pttabnew(n) |
---|
504 | nbtabold(n) = petabold(n) - pttabold(n) |
---|
505 | enddo |
---|
506 | ! |
---|
507 | eps = min(minval(dsnew(1:nbdim)),minval(dsold(1:nbdim))) / 100. |
---|
508 | ! |
---|
509 | do n = 1,nbdim |
---|
510 | if (snew(n) > (sold(n)+dsold(n)*nbtabold(n)+eps)) return |
---|
511 | if ((snew(n)+dsnew(n)*nbtabnew(n)-eps) < sold(n)) return |
---|
512 | enddo |
---|
513 | ! |
---|
514 | call Agrif_CopynD(new_var,old_var,pttabold,petabold,pttabnew,petabnew, & |
---|
515 | sold,snew,dsold,dsnew,nbdim) |
---|
516 | !--------------------------------------------------------------------------------------------------- |
---|
517 | end subroutine Agrif_Copy |
---|
518 | !=================================================================================================== |
---|
519 | ! |
---|
520 | !=================================================================================================== |
---|
521 | ! subroutine Agrif_CopynD |
---|
522 | ! |
---|
523 | !> Copy from the nD variable Old_Var the nD variable New_Var |
---|
524 | !--------------------------------------------------------------------------------------------------- |
---|
525 | subroutine Agrif_CopynD ( new_var, old_var, pttabold, petabold, pttabnew, petabnew, & |
---|
526 | sold, snew, dsold, dsnew, nbdim ) |
---|
527 | !--------------------------------------------------------------------------------------------------- |
---|
528 | type(Agrif_Variable), pointer, intent(inout) :: new_var |
---|
529 | type(Agrif_Variable), pointer, intent(in) :: old_var |
---|
530 | integer, dimension(nbdim), intent(in) :: pttabnew |
---|
531 | integer, dimension(nbdim), intent(in) :: petabnew |
---|
532 | integer, dimension(nbdim), intent(in) :: pttabold |
---|
533 | integer, dimension(nbdim), intent(in) :: petabold |
---|
534 | real, dimension(nbdim), intent(in) :: snew, sold |
---|
535 | real, dimension(nbdim), intent(in) :: dsnew,dsold |
---|
536 | integer, intent(in) :: nbdim |
---|
537 | ! |
---|
538 | integer :: i,j,k,l,m,n,i0,j0,k0,l0,m0,n0 |
---|
539 | ! |
---|
540 | real, dimension(nbdim) :: dim_gmin, dim_gmax |
---|
541 | real, dimension(nbdim) :: dim_newmin, dim_newmax |
---|
542 | real, dimension(nbdim) :: dim_min |
---|
543 | integer, dimension(nbdim) :: ind_gmin,ind_newmin, ind_newmax |
---|
544 | ! |
---|
545 | do i = 1,nbdim |
---|
546 | ! |
---|
547 | dim_gmin(i) = sold(i) |
---|
548 | dim_gmax(i) = sold(i) + (petabold(i)-pttabold(i)) * dsold(i) |
---|
549 | ! |
---|
550 | dim_newmin(i) = snew(i) |
---|
551 | dim_newmax(i) = snew(i) + (petabnew(i)-pttabnew(i)) * dsnew(i) |
---|
552 | ! |
---|
553 | enddo |
---|
554 | ! |
---|
555 | do i = 1,nbdim |
---|
556 | if (dim_gmax(i) < dim_newmin(i)) return |
---|
557 | if (dim_gmin(i) > dim_newmax(i)) return |
---|
558 | enddo |
---|
559 | ! |
---|
560 | do i = 1,nbdim |
---|
561 | ! |
---|
562 | ind_newmin(i) = pttabnew(i) - floor(-(max(dim_gmin(i),dim_newmin(i))-dim_newmin(i))/dsnew(i)) |
---|
563 | dim_min(i) = snew(i) + (ind_newmin(i)-pttabnew(i))*dsnew(i) |
---|
564 | ind_gmin(i) = pttabold(i) + nint((dim_min(i)-dim_gmin(i))/dsold(i)) |
---|
565 | ind_newmax(i) = pttabnew(i) + int((min(dim_gmax(i),dim_newmax(i))-dim_newmin(i))/dsnew(i)) |
---|
566 | ! |
---|
567 | enddo |
---|
568 | ! |
---|
569 | select case (nbdim) |
---|
570 | ! |
---|
571 | case (1) |
---|
572 | i0 = ind_gmin(1) |
---|
573 | do i = ind_newmin(1),ind_newmax(1) |
---|
574 | new_var % array1(i) = old_var % array1(i0) |
---|
575 | new_var % restore1D(i) = 1 |
---|
576 | i0 = i0 + int(dsnew(1)/dsold(1)) |
---|
577 | enddo |
---|
578 | ! |
---|
579 | case (2) |
---|
580 | i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1) |
---|
581 | j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2) |
---|
582 | new_var % array2(i,j) = old_var % array2(i0,j0) |
---|
583 | new_var % restore2D(i,j) = 1 |
---|
584 | j0 = j0 + int(dsnew(2)/dsold(2)) |
---|
585 | enddo |
---|
586 | i0 = i0 + int(dsnew(1)/dsold(1)) |
---|
587 | enddo |
---|
588 | ! |
---|
589 | case (3) |
---|
590 | i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1) |
---|
591 | j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2) |
---|
592 | k0 = ind_gmin(3) ; do k = ind_newmin(3),ind_newmax(3) |
---|
593 | new_var % array3(i,j,k) = old_var % array3(i0,j0,k0) |
---|
594 | new_var % restore3D(i,j,k) = 1 |
---|
595 | k0 = k0 + int(dsnew(3)/dsold(3)) |
---|
596 | enddo |
---|
597 | j0 = j0 + int(dsnew(2)/dsold(2)) |
---|
598 | enddo |
---|
599 | i0 = i0 + int(dsnew(1)/dsold(1)) |
---|
600 | enddo |
---|
601 | ! |
---|
602 | case (4) |
---|
603 | i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1) |
---|
604 | j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2) |
---|
605 | k0 = ind_gmin(3) ; do k = ind_newmin(3),ind_newmax(3) |
---|
606 | l0 = ind_gmin(4) ; do l = ind_newmin(4),ind_newmax(4) |
---|
607 | new_var % array4(i,j,k,l) = old_var % array4(i0,j0,k0,l0) |
---|
608 | new_var % restore4D(i,j,k,l) = 1 |
---|
609 | l0 = l0 + int(dsnew(4)/dsold(4)) |
---|
610 | enddo |
---|
611 | k0 = k0 + int(dsnew(3)/dsold(3)) |
---|
612 | enddo |
---|
613 | j0 = j0 + int(dsnew(2)/dsold(2)) |
---|
614 | enddo |
---|
615 | i0 = i0 + int(dsnew(1)/dsold(1)) |
---|
616 | enddo |
---|
617 | ! |
---|
618 | case (5) |
---|
619 | i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1) |
---|
620 | j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2) |
---|
621 | k0 = ind_gmin(3) ; do k = ind_newmin(3),ind_newmax(3) |
---|
622 | l0 = ind_gmin(4) ; do l = ind_newmin(4),ind_newmax(4) |
---|
623 | m0 = ind_gmin(5) ; do m = ind_newmin(5),ind_newmax(5) |
---|
624 | new_var % array5(i,j,k,l,m) = old_var % array5(i0,j0,k0,l0,m0) |
---|
625 | new_var % restore5D(i,j,k,l,m) = 1 |
---|
626 | m0 = m0 + int(dsnew(5)/dsold(5)) |
---|
627 | enddo |
---|
628 | l0 = l0 + int(dsnew(4)/dsold(4)) |
---|
629 | enddo |
---|
630 | k0 = k0 + int(dsnew(3)/dsold(3)) |
---|
631 | enddo |
---|
632 | j0 = j0 + int(dsnew(2)/dsold(2)) |
---|
633 | enddo |
---|
634 | i0 = i0 + int(dsnew(1)/dsold(1)) |
---|
635 | enddo |
---|
636 | ! |
---|
637 | case (6) |
---|
638 | i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1) |
---|
639 | j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2) |
---|
640 | k0 = ind_gmin(3) ; do k = ind_newmin(3),ind_newmax(3) |
---|
641 | l0 = ind_gmin(4) ; do l = ind_newmin(4),ind_newmax(4) |
---|
642 | m0 = ind_gmin(5) ; do m = ind_newmin(5),ind_newmax(5) |
---|
643 | n0 = ind_gmin(6) ; do n = ind_newmin(6),ind_newmax(6) |
---|
644 | new_var % array6(i,j,k,l,m,n) = old_var % array6(i0,j0,k0,l0,m0,n0) |
---|
645 | new_var % restore6D(i,j,k,l,m,n) = 1 |
---|
646 | n0 = n0 + int(dsnew(6)/dsold(6)) |
---|
647 | enddo |
---|
648 | m0 = m0 + int(dsnew(5)/dsold(5)) |
---|
649 | enddo |
---|
650 | l0 = l0 + int(dsnew(4)/dsold(4)) |
---|
651 | enddo |
---|
652 | k0 = k0 + int(dsnew(3)/dsold(3)) |
---|
653 | enddo |
---|
654 | j0 = j0 + int(dsnew(2)/dsold(2)) |
---|
655 | enddo |
---|
656 | i0 = i0 + int(dsnew(1)/dsold(1)) |
---|
657 | enddo |
---|
658 | ! |
---|
659 | end select |
---|
660 | !--------------------------------------------------------------------------------------------------- |
---|
661 | end subroutine Agrif_CopynD |
---|
662 | !=================================================================================================== |
---|
663 | ! |
---|
664 | end module Agrif_Save |
---|