1 | ! |
---|
2 | ! $Id$ |
---|
3 | ! |
---|
4 | ! AGRIF (Adaptive Grid Refinement In Fortran) |
---|
5 | ! |
---|
6 | ! Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr) |
---|
7 | ! Christophe Vouland (Christophe.Vouland@imag.fr) |
---|
8 | ! |
---|
9 | ! This program is free software; you can redistribute it and/or modify |
---|
10 | ! it under the terms of the GNU General Public License as published by |
---|
11 | ! the Free Software Foundation; either version 2 of the License, or |
---|
12 | ! (at your option) any later version. |
---|
13 | ! |
---|
14 | ! This program is distributed in the hope that it will be useful, |
---|
15 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
16 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
17 | ! GNU General Public License for more details. |
---|
18 | ! |
---|
19 | ! You should have received a copy of the GNU General Public License |
---|
20 | ! along with this program; if not, write to the Free Software |
---|
21 | ! Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
---|
22 | ! |
---|
23 | !> Module Agrif_BcFunction. |
---|
24 | ! |
---|
25 | module Agrif_BcFunction |
---|
26 | ! |
---|
27 | ! Modules used: |
---|
28 | ! |
---|
29 | use Agrif_Boundary |
---|
30 | use Agrif_Update |
---|
31 | use Agrif_Save |
---|
32 | ! |
---|
33 | implicit none |
---|
34 | ! |
---|
35 | interface Agrif_Set_Parent |
---|
36 | module procedure Agrif_Set_Parent_int, & |
---|
37 | Agrif_Set_Parent_real4, & |
---|
38 | Agrif_Set_Parent_real8 |
---|
39 | end interface |
---|
40 | ! |
---|
41 | interface Agrif_Save_Forrestore |
---|
42 | module procedure Agrif_Save_Forrestore0d, & |
---|
43 | Agrif_Save_Forrestore2d, & |
---|
44 | Agrif_Save_Forrestore3d, & |
---|
45 | Agrif_Save_Forrestore4d |
---|
46 | end interface |
---|
47 | ! |
---|
48 | contains |
---|
49 | ! |
---|
50 | !=================================================================================================== |
---|
51 | ! subroutine Agrif_Set_parent_int |
---|
52 | ! |
---|
53 | !> To set the TYPE of the variable |
---|
54 | !--------------------------------------------------------------------------------------------------- |
---|
55 | subroutine Agrif_Set_parent_int(integer_variable,value) |
---|
56 | !--------------------------------------------------------------------------------------------------- |
---|
57 | integer, intent(in) :: integer_variable !< indice of the variable in tabvars |
---|
58 | integer, intent(in) :: value !< input value |
---|
59 | ! |
---|
60 | |
---|
61 | integer :: i |
---|
62 | logical :: i_found |
---|
63 | |
---|
64 | i_found = .FALSE. |
---|
65 | |
---|
66 | do i=1,Agrif_NbVariables(4) |
---|
67 | if (LOC(integer_variable) == LOC(agrif_curgrid%tabvars_i(i)%iarray0)) then |
---|
68 | agrif_curgrid%tabvars_i(i)%parent_var%iarray0 = value |
---|
69 | i_found = .TRUE. |
---|
70 | EXIT |
---|
71 | endif |
---|
72 | enddo |
---|
73 | |
---|
74 | if (.NOT.i_found) STOP 'Agrif_Set_Integer : Variable not found' |
---|
75 | |
---|
76 | !--------------------------------------------------------------------------------------------------- |
---|
77 | end subroutine Agrif_Set_parent_int |
---|
78 | !=================================================================================================== |
---|
79 | ! |
---|
80 | !=================================================================================================== |
---|
81 | ! subroutine Agrif_Set_parent_real4 |
---|
82 | !--------------------------------------------------------------------------------------------------- |
---|
83 | !> To set the parent value of a real variable |
---|
84 | !--------------------------------------------------------------------------------------------------- |
---|
85 | subroutine Agrif_Set_parent_real4 ( real_variable, value ) |
---|
86 | !--------------------------------------------------------------------------------------------------- |
---|
87 | real(kind=4), intent(in) :: real_variable !< input variable |
---|
88 | real(kind=4),intent(in) :: value !< input value for the parent grid |
---|
89 | |
---|
90 | integer :: i |
---|
91 | logical :: i_found |
---|
92 | |
---|
93 | i_found = .FALSE. |
---|
94 | |
---|
95 | do i=1,Agrif_NbVariables(2) |
---|
96 | if (LOC(real_variable) == LOC(agrif_curgrid%tabvars_r(i)%array0)) then |
---|
97 | agrif_curgrid%tabvars_r(i)%parent_var%array0 = value |
---|
98 | agrif_curgrid%tabvars_r(i)%parent_var%sarray0 = value |
---|
99 | i_found = .TRUE. |
---|
100 | EXIT |
---|
101 | endif |
---|
102 | enddo |
---|
103 | |
---|
104 | IF (.NOT.i_found) THEN |
---|
105 | do i=1,Agrif_NbVariables(2) |
---|
106 | if (LOC(real_variable) == LOC(agrif_curgrid%tabvars_r(i)%sarray0)) then |
---|
107 | agrif_curgrid%tabvars_r(i)%parent_var%array0 = value |
---|
108 | agrif_curgrid%tabvars_r(i)%parent_var%sarray0 = value |
---|
109 | i_found = .TRUE. |
---|
110 | EXIT |
---|
111 | endif |
---|
112 | enddo |
---|
113 | ENDIF |
---|
114 | |
---|
115 | if (.NOT.i_found) STOP 'Agrif_Set_parent_real4 : Variable not found' |
---|
116 | |
---|
117 | !--------------------------------------------------------------------------------------------------- |
---|
118 | end subroutine Agrif_Set_parent_real4 |
---|
119 | !=================================================================================================== |
---|
120 | ! |
---|
121 | !=================================================================================================== |
---|
122 | ! subroutine Agrif_Set_parent_real8 |
---|
123 | !--------------------------------------------------------------------------------------------------- |
---|
124 | !> To set the parent value of a real variable |
---|
125 | !--------------------------------------------------------------------------------------------------- |
---|
126 | subroutine Agrif_Set_parent_real8 ( real_variable, value ) |
---|
127 | !--------------------------------------------------------------------------------------------------- |
---|
128 | real(kind=8), intent(in) :: real_variable !< input variable |
---|
129 | real(kind=8),intent(in) :: value !< input value for the parent grid |
---|
130 | |
---|
131 | integer :: i |
---|
132 | logical :: i_found |
---|
133 | |
---|
134 | i_found = .FALSE. |
---|
135 | |
---|
136 | do i=1,Agrif_NbVariables(2) |
---|
137 | if (LOC(real_variable) == LOC(agrif_curgrid%tabvars_r(i)%array0)) then |
---|
138 | agrif_curgrid%tabvars_r(i)%parent_var%darray0 = value |
---|
139 | agrif_curgrid%tabvars_r(i)%parent_var%array0 = value |
---|
140 | i_found = .TRUE. |
---|
141 | EXIT |
---|
142 | endif |
---|
143 | enddo |
---|
144 | |
---|
145 | IF (.NOT.i_found) THEN |
---|
146 | do i=1,Agrif_NbVariables(2) |
---|
147 | if (LOC(real_variable) == LOC(agrif_curgrid%tabvars_r(i)%darray0)) then |
---|
148 | agrif_curgrid%tabvars_r(i)%parent_var%darray0 = value |
---|
149 | agrif_curgrid%tabvars_r(i)%parent_var%array0 = value |
---|
150 | i_found = .TRUE. |
---|
151 | EXIT |
---|
152 | endif |
---|
153 | enddo |
---|
154 | ENDIF |
---|
155 | |
---|
156 | if (.NOT.i_found) STOP 'Agrif_Set_parent_real8 : Variable not found' |
---|
157 | |
---|
158 | !--------------------------------------------------------------------------------------------------- |
---|
159 | end subroutine Agrif_Set_parent_real8 |
---|
160 | !=================================================================================================== |
---|
161 | ! |
---|
162 | !=================================================================================================== |
---|
163 | ! subroutine Agrif_Set_bc |
---|
164 | !--------------------------------------------------------------------------------------------------- |
---|
165 | subroutine Agrif_Set_bc ( tabvarsindic, bcinfsup, Interpolationshouldbemade ) |
---|
166 | !--------------------------------------------------------------------------------------------------- |
---|
167 | integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars |
---|
168 | integer, dimension(2), intent(in) :: bcinfsup !< bcinfsup |
---|
169 | logical, optional, intent(in) :: Interpolationshouldbemade !< interpolation should be made |
---|
170 | ! |
---|
171 | integer :: indic ! indice of the variable in tabvars |
---|
172 | type(Agrif_Variable), pointer :: var |
---|
173 | ! |
---|
174 | var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
175 | if (.not.associated(var)) return ! Grand mother grid case |
---|
176 | ! |
---|
177 | if ( Agrif_Curgrid % fixedrank /= 0 ) then |
---|
178 | if ( .not.associated(var % oldvalues2D) ) then |
---|
179 | allocate(var % oldvalues2D(2,1)) |
---|
180 | var % interpIndex = -1 |
---|
181 | var % oldvalues2D = 0. |
---|
182 | endif |
---|
183 | if ( present(Interpolationshouldbemade) ) then |
---|
184 | var % Interpolationshouldbemade = Interpolationshouldbemade |
---|
185 | endif |
---|
186 | endif |
---|
187 | ! |
---|
188 | var % bcinf = bcinfsup(1) |
---|
189 | var % bcsup = bcinfsup(2) |
---|
190 | !--------------------------------------------------------------------------------------------------- |
---|
191 | end subroutine Agrif_Set_bc |
---|
192 | !=================================================================================================== |
---|
193 | ! |
---|
194 | !=================================================================================================== |
---|
195 | ! subroutine Agrif_Set_interp |
---|
196 | !--------------------------------------------------------------------------------------------------- |
---|
197 | subroutine Agrif_Set_interp ( tabvarsindic, interp, interp1, interp2, interp3 , interp4) |
---|
198 | !--------------------------------------------------------------------------------------------------- |
---|
199 | integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars |
---|
200 | integer, optional, intent(in) :: interp, interp1, interp2, interp3, interp4 |
---|
201 | ! |
---|
202 | integer :: indic ! indice of the variable in tabvars |
---|
203 | type(Agrif_Variable), pointer :: var |
---|
204 | ! |
---|
205 | var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
206 | if (.not.associated(var)) return ! Grand mother grid case |
---|
207 | ! |
---|
208 | var % type_interp = Agrif_Constant |
---|
209 | ! |
---|
210 | if (present(interp)) var % type_interp = interp |
---|
211 | if (present(interp1)) var % type_interp(1) = interp1 |
---|
212 | if (present(interp2)) var % type_interp(2) = interp2 |
---|
213 | if (present(interp3)) var % type_interp(3) = interp3 |
---|
214 | if (present(interp4)) var % type_interp(4) = interp4 |
---|
215 | !--------------------------------------------------------------------------------------------------- |
---|
216 | end subroutine Agrif_Set_interp |
---|
217 | !=================================================================================================== |
---|
218 | ! |
---|
219 | !=================================================================================================== |
---|
220 | ! subroutine Agrif_Set_bcinterp |
---|
221 | !--------------------------------------------------------------------------------------------------- |
---|
222 | subroutine Agrif_Set_bcinterp ( tabvarsindic, interp, interp1, interp2, interp3, interp4, & |
---|
223 | interp11, interp12, interp21, interp22 ) |
---|
224 | !--------------------------------------------------------------------------------------------------- |
---|
225 | INTEGER, intent(in) :: tabvarsindic !< indice of the variable in tabvars |
---|
226 | INTEGER, OPTIONAL, intent(in) :: interp, interp1, interp2, interp3, interp4 |
---|
227 | INTEGER, OPTIONAL, intent(in) :: interp11, interp12, interp21, interp22 |
---|
228 | ! |
---|
229 | INTEGER :: indic ! indice of the variable in tabvars |
---|
230 | TYPE(Agrif_Variable), pointer :: var |
---|
231 | ! |
---|
232 | var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
233 | ! |
---|
234 | var % type_interp_bc = Agrif_Constant |
---|
235 | ! |
---|
236 | if (present(interp)) var % type_interp_bc = interp |
---|
237 | if (present(interp1)) var % type_interp_bc(:,1) = interp1 |
---|
238 | if (present(interp11)) var % type_interp_bc(1,1) = interp11 |
---|
239 | if (present(interp12)) var % type_interp_bc(1,2) = interp12 |
---|
240 | if (present(interp2)) var % type_interp_bc(:,2) = interp2 |
---|
241 | if (present(interp21)) var % type_interp_bc(2,1) = interp21 |
---|
242 | if (present(interp22)) var % type_interp_bc(2,2) = interp22 |
---|
243 | if (present(interp3)) var % type_interp_bc(:,3) = interp3 |
---|
244 | if (present(interp4)) var % type_interp_bc(:,4) = interp4 |
---|
245 | !--------------------------------------------------------------------------------------------------- |
---|
246 | end subroutine Agrif_Set_bcinterp |
---|
247 | !=================================================================================================== |
---|
248 | ! |
---|
249 | !=================================================================================================== |
---|
250 | ! subroutine Agrif_Set_UpdateType |
---|
251 | !--------------------------------------------------------------------------------------------------- |
---|
252 | subroutine Agrif_Set_UpdateType ( tabvarsindic, update, update1, update2, & |
---|
253 | update3, update4, update5 ) |
---|
254 | !--------------------------------------------------------------------------------------------------- |
---|
255 | INTEGER, intent(in) :: tabvarsindic !< indice of the variable in tabvars |
---|
256 | INTEGER, OPTIONAL, intent(in) :: update, update1, update2, update3, update4, update5 |
---|
257 | ! |
---|
258 | INTEGER :: indic ! indice of the variable in tabvars |
---|
259 | type(Agrif_Variable), pointer :: root_var |
---|
260 | ! |
---|
261 | |
---|
262 | root_var => Agrif_Search_Variable(Agrif_Mygrid,tabvarsindic) |
---|
263 | |
---|
264 | ! |
---|
265 | root_var % type_update = Agrif_Update_Copy |
---|
266 | if (present(update)) root_var % type_update = update |
---|
267 | if (present(update1)) root_var % type_update(1) = update1 |
---|
268 | if (present(update2)) root_var % type_update(2) = update2 |
---|
269 | if (present(update3)) root_var % type_update(3) = update3 |
---|
270 | if (present(update4)) root_var % type_update(4) = update4 |
---|
271 | if (present(update5)) root_var % type_update(5) = update5 |
---|
272 | !--------------------------------------------------------------------------------------------------- |
---|
273 | end subroutine Agrif_Set_UpdateType |
---|
274 | !=================================================================================================== |
---|
275 | ! |
---|
276 | !=================================================================================================== |
---|
277 | ! subroutine Agrif_Set_restore |
---|
278 | !--------------------------------------------------------------------------------------------------- |
---|
279 | subroutine Agrif_Set_restore ( tabvarsindic ) |
---|
280 | !--------------------------------------------------------------------------------------------------- |
---|
281 | INTEGER, intent(in) :: tabvarsindic !< indice of the variable in tabvars |
---|
282 | ! |
---|
283 | INTEGER :: indic ! indice of the variable in tabvars |
---|
284 | ! |
---|
285 | print *,'CURRENTLY BROKEN' |
---|
286 | STOP |
---|
287 | |
---|
288 | indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0 |
---|
289 | ! |
---|
290 | Agrif_Mygrid%tabvars(indic) % restore = .TRUE. |
---|
291 | !--------------------------------------------------------------------------------------------------- |
---|
292 | end subroutine Agrif_Set_restore |
---|
293 | !=================================================================================================== |
---|
294 | ! |
---|
295 | !=================================================================================================== |
---|
296 | ! subroutine Agrif_Init_variable |
---|
297 | !--------------------------------------------------------------------------------------------------- |
---|
298 | subroutine Agrif_Init_variable ( tabvarsindic, procname ) |
---|
299 | !--------------------------------------------------------------------------------------------------- |
---|
300 | INTEGER, intent(in) :: tabvarsindic !< indice of the variable in tabvars |
---|
301 | procedure() :: procname !< Data recovery procedure |
---|
302 | ! |
---|
303 | if ( Agrif_Curgrid%level <= 0 ) return |
---|
304 | ! |
---|
305 | call Agrif_Interp_variable(tabvarsindic, procname) |
---|
306 | call Agrif_Bc_variable(tabvarsindic, procname, 1.) |
---|
307 | !--------------------------------------------------------------------------------------------------- |
---|
308 | end subroutine Agrif_Init_variable |
---|
309 | !=================================================================================================== |
---|
310 | ! |
---|
311 | !=================================================================================================== |
---|
312 | ! subroutine Agrif_Bc_variable |
---|
313 | !--------------------------------------------------------------------------------------------------- |
---|
314 | subroutine Agrif_Bc_variable ( tabvarsindic, procname, calledweight ) |
---|
315 | !--------------------------------------------------------------------------------------------------- |
---|
316 | integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars |
---|
317 | procedure() :: procname |
---|
318 | real, optional, intent(in) :: calledweight |
---|
319 | ! |
---|
320 | real :: weight |
---|
321 | logical :: pweight |
---|
322 | integer :: indic |
---|
323 | integer :: nbdim |
---|
324 | type(Agrif_Variable), pointer :: root_var |
---|
325 | type(Agrif_Variable), pointer :: parent_var |
---|
326 | type(Agrif_Variable), pointer :: child_var |
---|
327 | type(Agrif_Variable), pointer :: child_tmp ! Temporary variable on the child grid |
---|
328 | integer :: i |
---|
329 | integer,dimension(7) :: lb, ub |
---|
330 | ! |
---|
331 | if ( Agrif_Curgrid%level <= 0 ) return |
---|
332 | ! |
---|
333 | ! |
---|
334 | if ( present(calledweight) ) then |
---|
335 | weight = calledweight |
---|
336 | pweight = .true. |
---|
337 | else |
---|
338 | weight = 0. |
---|
339 | pweight = .false. |
---|
340 | endif |
---|
341 | ! |
---|
342 | child_var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
343 | parent_var => child_var % parent_var |
---|
344 | root_var => child_var % root_var |
---|
345 | ! |
---|
346 | nbdim = root_var % nbdim |
---|
347 | ! |
---|
348 | do i=1,nbdim |
---|
349 | if (root_var%coords(i) == 0) then |
---|
350 | lb(i) = parent_var%lb(i) |
---|
351 | ub(i) = parent_var%ub(i) |
---|
352 | else |
---|
353 | lb(i) = child_var%lb(i) |
---|
354 | ub(i) = child_var%ub(i) |
---|
355 | endif |
---|
356 | enddo |
---|
357 | |
---|
358 | select case( nbdim ) |
---|
359 | case(1) |
---|
360 | allocate(parray1(lb(1):ub(1))) |
---|
361 | case(2) |
---|
362 | allocate(parray2(lb(1):ub(1), & |
---|
363 | lb(2):ub(2) )) |
---|
364 | case(3) |
---|
365 | allocate(parray3(lb(1):ub(1), & |
---|
366 | lb(2):ub(2), & |
---|
367 | lb(3):ub(3) )) |
---|
368 | case(4) |
---|
369 | allocate(parray4(lb(1):ub(1), & |
---|
370 | lb(2):ub(2), & |
---|
371 | lb(3):ub(3), & |
---|
372 | lb(4):ub(4) )) |
---|
373 | case(5) |
---|
374 | allocate(parray5(lb(1):ub(1), & |
---|
375 | lb(2):ub(2), & |
---|
376 | lb(3):ub(3), & |
---|
377 | lb(4):ub(4), & |
---|
378 | lb(5):ub(5) )) |
---|
379 | case(6) |
---|
380 | allocate(parray6(lb(1):ub(1), & |
---|
381 | lb(2):ub(2), & |
---|
382 | lb(3):ub(3), & |
---|
383 | lb(4):ub(4), & |
---|
384 | lb(5):ub(5), & |
---|
385 | lb(6):ub(6) )) |
---|
386 | end select |
---|
387 | ! |
---|
388 | ! Create temporary child variable |
---|
389 | allocate(child_tmp) |
---|
390 | ! |
---|
391 | child_tmp % root_var => root_var |
---|
392 | child_tmp % parent_var => parent_var |
---|
393 | child_tmp % oldvalues2D => child_var % oldvalues2D |
---|
394 | ! |
---|
395 | ! Index indicating if a space interpolation is necessary |
---|
396 | child_tmp % interpIndex = child_var % interpIndex |
---|
397 | child_tmp % list_interp => child_var % list_interp |
---|
398 | child_tmp % Interpolationshouldbemade = child_var % Interpolationshouldbemade |
---|
399 | ! |
---|
400 | child_tmp % point = child_var % point |
---|
401 | child_tmp % lb = child_var % lb |
---|
402 | child_tmp % ub = child_var % ub |
---|
403 | ! |
---|
404 | child_tmp % bcinf = child_var % bcinf |
---|
405 | child_tmp % bcsup = child_var % bcsup |
---|
406 | ! |
---|
407 | child_tmp % childarray = child_var % childarray |
---|
408 | child_tmp % memberin = child_var % memberin |
---|
409 | ! |
---|
410 | call Agrif_CorrectVariable(parent_var, child_tmp, pweight, weight, procname) |
---|
411 | ! |
---|
412 | child_var % childarray = child_tmp % childarray |
---|
413 | child_var % memberin = child_tmp % memberin |
---|
414 | ! |
---|
415 | child_var % oldvalues2D => child_tmp % oldvalues2D |
---|
416 | child_var % list_interp => child_tmp % list_interp |
---|
417 | ! |
---|
418 | child_var % interpIndex = child_tmp % interpIndex |
---|
419 | ! |
---|
420 | deallocate(child_tmp) |
---|
421 | ! |
---|
422 | select case( nbdim ) |
---|
423 | case(1); deallocate(parray1) |
---|
424 | case(2); deallocate(parray2) |
---|
425 | case(3); deallocate(parray3) |
---|
426 | case(4); deallocate(parray4) |
---|
427 | case(5); deallocate(parray5) |
---|
428 | case(6); deallocate(parray6) |
---|
429 | end select |
---|
430 | !--------------------------------------------------------------------------------------------------- |
---|
431 | end subroutine Agrif_Bc_variable |
---|
432 | !=================================================================================================== |
---|
433 | ! |
---|
434 | !=================================================================================================== |
---|
435 | ! subroutine Agrif_Interp_variable |
---|
436 | !--------------------------------------------------------------------------------------------------- |
---|
437 | subroutine Agrif_Interp_variable ( tabvarsindic, procname ) |
---|
438 | !--------------------------------------------------------------------------------------------------- |
---|
439 | integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars |
---|
440 | procedure() :: procname !< Data recovery procedure |
---|
441 | ! |
---|
442 | integer :: nbdim |
---|
443 | integer :: indic ! indice of the variable in tabvars |
---|
444 | logical :: torestore |
---|
445 | type(Agrif_Variable), pointer :: root_var |
---|
446 | type(Agrif_Variable), pointer :: parent_var ! Variable on the parent grid |
---|
447 | type(Agrif_Variable), pointer :: child_var ! Variable on the parent grid |
---|
448 | type(Agrif_Variable), pointer :: child_tmp ! Temporary variable on the child grid |
---|
449 | ! |
---|
450 | |
---|
451 | if ( Agrif_Curgrid%level <= 0 ) return |
---|
452 | ! |
---|
453 | |
---|
454 | child_var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
455 | parent_var => child_var % parent_var |
---|
456 | root_var => child_var % root_var |
---|
457 | |
---|
458 | ! |
---|
459 | nbdim = root_var % nbdim |
---|
460 | torestore = root_var % restore |
---|
461 | ! |
---|
462 | allocate(child_tmp) |
---|
463 | ! |
---|
464 | child_tmp % root_var => root_var |
---|
465 | child_tmp % parent_var => parent_var |
---|
466 | child_tmp % nbdim = root_var % nbdim |
---|
467 | child_tmp % point = child_var % point |
---|
468 | child_tmp % lb = child_var % lb |
---|
469 | child_tmp % ub = child_var % ub |
---|
470 | child_tmp % interpIndex = child_var % interpIndex |
---|
471 | child_tmp % list_interp => child_var % list_interp |
---|
472 | child_tmp % Interpolationshouldbemade = child_var % Interpolationshouldbemade |
---|
473 | ! |
---|
474 | if ( torestore ) then |
---|
475 | select case( nbdim ) |
---|
476 | case(1) |
---|
477 | parray1 = child_var % array1 |
---|
478 | child_tmp % restore1D => child_var % restore1D |
---|
479 | case(2) |
---|
480 | parray2 = child_var % array2 |
---|
481 | child_tmp % restore2D => child_var % restore2D |
---|
482 | case(3) |
---|
483 | parray3 = child_var % array3 |
---|
484 | child_tmp % restore3D => child_var % restore3D |
---|
485 | case(4) |
---|
486 | parray4 = child_var % array4 |
---|
487 | child_tmp % restore4D => child_var % restore4D |
---|
488 | case(5) |
---|
489 | parray5 = child_var % array5 |
---|
490 | child_tmp % restore5D => child_var % restore5D |
---|
491 | case(6) |
---|
492 | parray6 = child_var % array6 |
---|
493 | child_tmp % restore6D => child_var % restore6D |
---|
494 | end select |
---|
495 | endif |
---|
496 | ! |
---|
497 | call Agrif_InterpVariable(parent_var, child_tmp, torestore, procname) |
---|
498 | ! |
---|
499 | child_var % list_interp => child_tmp % list_interp |
---|
500 | ! |
---|
501 | deallocate(child_tmp) |
---|
502 | !--------------------------------------------------------------------------------------------------- |
---|
503 | end subroutine Agrif_Interp_variable |
---|
504 | !=================================================================================================== |
---|
505 | ! |
---|
506 | !=================================================================================================== |
---|
507 | ! subroutine Agrif_Update_Variable |
---|
508 | !--------------------------------------------------------------------------------------------------- |
---|
509 | subroutine Agrif_Update_Variable ( tabvarsindic, procname, & |
---|
510 | locupdate, locupdate1, locupdate2, locupdate3, locupdate4 ) |
---|
511 | !--------------------------------------------------------------------------------------------------- |
---|
512 | integer, intent(in) :: tabvarsindic !< Indice of the variable in tabvars |
---|
513 | procedure() :: procname !< Data recovery procedure |
---|
514 | integer, dimension(2), intent(in), optional :: locupdate |
---|
515 | integer, dimension(2), intent(in), optional :: locupdate1 |
---|
516 | integer, dimension(2), intent(in), optional :: locupdate2 |
---|
517 | integer, dimension(2), intent(in), optional :: locupdate3 |
---|
518 | integer, dimension(2), intent(in), optional :: locupdate4 |
---|
519 | !--------------------------------------------------------------------------------------------------- |
---|
520 | integer :: indic |
---|
521 | integer :: nbdim |
---|
522 | integer, dimension(6) :: updateinf ! First positions where interpolations are calculated |
---|
523 | integer, dimension(6) :: updatesup ! Last positions where interpolations are calculated |
---|
524 | type(Agrif_Variable), pointer :: root_var |
---|
525 | type(Agrif_Variable), pointer :: parent_var |
---|
526 | type(Agrif_Variable), pointer :: child_var |
---|
527 | ! |
---|
528 | if ( Agrif_Root() .AND. (.not.agrif_coarse) ) return |
---|
529 | if (agrif_curgrid%grand_mother_grid) return |
---|
530 | ! |
---|
531 | |
---|
532 | child_var => Agrif_Search_Variable(Agrif_Curgrid, tabvarsindic) |
---|
533 | parent_var => child_var % parent_var |
---|
534 | |
---|
535 | if (.not.associated(parent_var)) then |
---|
536 | ! can occur during the first update of Agrif_Coarsegrid (if any) |
---|
537 | parent_var => Agrif_Search_Variable(Agrif_Curgrid % parent, tabvarsindic) |
---|
538 | child_var % parent_var => parent_var |
---|
539 | endif |
---|
540 | |
---|
541 | root_var => child_var % root_var |
---|
542 | |
---|
543 | ! |
---|
544 | nbdim = root_var % nbdim |
---|
545 | ! |
---|
546 | updateinf = -99 |
---|
547 | updatesup = -99 |
---|
548 | ! |
---|
549 | if ( present(locupdate) ) then |
---|
550 | updateinf(1:nbdim) = locupdate(1) |
---|
551 | updatesup(1:nbdim) = locupdate(2) |
---|
552 | endif |
---|
553 | ! |
---|
554 | if ( present(locupdate1) ) then |
---|
555 | updateinf(1) = locupdate1(1) |
---|
556 | updatesup(1) = locupdate1(2) |
---|
557 | endif |
---|
558 | ! |
---|
559 | if ( present(locupdate2) ) then |
---|
560 | updateinf(2) = locupdate2(1) |
---|
561 | updatesup(2) = locupdate2(2) |
---|
562 | endif |
---|
563 | |
---|
564 | if ( present(locupdate3) ) then |
---|
565 | updateinf(3) = locupdate3(1) |
---|
566 | updatesup(3) = locupdate3(2) |
---|
567 | endif |
---|
568 | |
---|
569 | if ( present(locupdate4) ) then |
---|
570 | updateinf(4) = locupdate4(1) |
---|
571 | updatesup(4) = locupdate4(2) |
---|
572 | endif |
---|
573 | ! |
---|
574 | call Agrif_UpdateVariable( parent_var, child_var, updateinf, updatesup, procname ) |
---|
575 | !--------------------------------------------------------------------------------------------------- |
---|
576 | end subroutine Agrif_Update_Variable |
---|
577 | !=================================================================================================== |
---|
578 | ! |
---|
579 | !=================================================================================================== |
---|
580 | ! subroutine Agrif_Save_ForRestore0D |
---|
581 | !--------------------------------------------------------------------------------------------------- |
---|
582 | subroutine Agrif_Save_ForRestore0D ( tabvarsindic0, tabvarsindic ) |
---|
583 | !--------------------------------------------------------------------------------------------------- |
---|
584 | integer, intent(in) :: tabvarsindic0, tabvarsindic |
---|
585 | ! |
---|
586 | type(Agrif_Variable), pointer :: root_var, save_var |
---|
587 | integer :: nbdim |
---|
588 | ! |
---|
589 | print *,'CURRENTLY BROKEN' |
---|
590 | STOP |
---|
591 | root_var => Agrif_Mygrid % tabvars(tabvarsindic0) |
---|
592 | save_var => Agrif_Curgrid % tabvars(tabvarsindic0) |
---|
593 | nbdim = root_var % nbdim |
---|
594 | ! |
---|
595 | select case(nbdim) |
---|
596 | case(2); call Agrif_Save_ForRestore2D(save_var % array2, tabvarsindic) |
---|
597 | case(3); call Agrif_Save_ForRestore3D(save_var % array3, tabvarsindic) |
---|
598 | case(4); call Agrif_Save_ForRestore4D(save_var % array4, tabvarsindic) |
---|
599 | end select |
---|
600 | !--------------------------------------------------------------------------------------------------- |
---|
601 | end subroutine Agrif_Save_ForRestore0D |
---|
602 | !=================================================================================================== |
---|
603 | ! |
---|
604 | !=================================================================================================== |
---|
605 | ! subroutine Agrif_Save_ForRestore2D |
---|
606 | !--------------------------------------------------------------------------------------------------- |
---|
607 | subroutine Agrif_Save_ForRestore2D ( q, tabvarsindic ) |
---|
608 | !--------------------------------------------------------------------------------------------------- |
---|
609 | real, dimension(:,:), intent(in) :: q |
---|
610 | integer, intent(in) :: tabvarsindic |
---|
611 | ! |
---|
612 | type(Agrif_Variable), pointer :: root_var, save_var |
---|
613 | integer :: indic |
---|
614 | ! |
---|
615 | print *,'CURRENTLY BROKEN' |
---|
616 | STOP |
---|
617 | indic = tabvarsindic |
---|
618 | if (tabvarsindic >= 0) then |
---|
619 | if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then |
---|
620 | indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0 |
---|
621 | endif |
---|
622 | endif |
---|
623 | ! |
---|
624 | if (indic <= 0) then |
---|
625 | save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic) |
---|
626 | root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic) |
---|
627 | else |
---|
628 | save_var => Agrif_Curgrid % tabvars(indic) |
---|
629 | root_var => Agrif_Mygrid % tabvars(indic) |
---|
630 | endif |
---|
631 | ! |
---|
632 | if ( .not.allocated(save_var%array2) ) then |
---|
633 | allocate(save_var%array2(save_var%lb(1):save_var%ub(1), & |
---|
634 | save_var%lb(2):save_var%ub(2))) |
---|
635 | endif |
---|
636 | ! |
---|
637 | save_var % array2 = q |
---|
638 | root_var % restore = .true. |
---|
639 | !--------------------------------------------------------------------------------------------------- |
---|
640 | end subroutine Agrif_Save_ForRestore2D |
---|
641 | !=================================================================================================== |
---|
642 | ! |
---|
643 | !=================================================================================================== |
---|
644 | ! subroutine Agrif_Save_ForRestore3D |
---|
645 | !--------------------------------------------------------------------------------------------------- |
---|
646 | subroutine Agrif_Save_ForRestore3D ( q, tabvarsindic ) |
---|
647 | !--------------------------------------------------------------------------------------------------- |
---|
648 | real, dimension(:,:,:), intent(in) :: q |
---|
649 | integer, intent(in) :: tabvarsindic |
---|
650 | ! |
---|
651 | type(Agrif_Variable), pointer :: root_var, save_var |
---|
652 | integer :: indic |
---|
653 | ! |
---|
654 | print *,'CURRENTLY BROKEN' |
---|
655 | STOP |
---|
656 | |
---|
657 | indic = tabvarsindic |
---|
658 | if (tabvarsindic >= 0) then |
---|
659 | if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then |
---|
660 | indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0 |
---|
661 | endif |
---|
662 | endif |
---|
663 | ! |
---|
664 | if (indic <= 0) then |
---|
665 | save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic) |
---|
666 | root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic) |
---|
667 | else |
---|
668 | save_var => Agrif_Curgrid % tabvars(indic) |
---|
669 | root_var => Agrif_Mygrid % tabvars(indic) |
---|
670 | endif |
---|
671 | ! |
---|
672 | if ( .not.allocated(save_var%array3) ) then |
---|
673 | allocate(save_var%array3(save_var%lb(1):save_var%ub(1), & |
---|
674 | save_var%lb(2):save_var%ub(2), & |
---|
675 | save_var%lb(3):save_var%ub(3))) |
---|
676 | endif |
---|
677 | ! |
---|
678 | save_var % array3 = q |
---|
679 | root_var % restore = .true. |
---|
680 | !--------------------------------------------------------------------------------------------------- |
---|
681 | end subroutine Agrif_Save_ForRestore3D |
---|
682 | !=================================================================================================== |
---|
683 | ! |
---|
684 | !=================================================================================================== |
---|
685 | ! subroutine Agrif_Save_ForRestore4D |
---|
686 | !--------------------------------------------------------------------------------------------------- |
---|
687 | subroutine Agrif_Save_ForRestore4D ( q, tabvarsindic ) |
---|
688 | !--------------------------------------------------------------------------------------------------- |
---|
689 | real, dimension(:,:,:,:), intent(in) :: q |
---|
690 | integer, intent(in) :: tabvarsindic |
---|
691 | ! |
---|
692 | type(Agrif_Variable), pointer :: root_var, save_var |
---|
693 | integer :: indic |
---|
694 | ! |
---|
695 | print *,'CURRENTLY BROKEN' |
---|
696 | STOP |
---|
697 | indic = tabvarsindic |
---|
698 | if (tabvarsindic >= 0) then |
---|
699 | if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then |
---|
700 | indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0 |
---|
701 | endif |
---|
702 | endif |
---|
703 | ! |
---|
704 | if (indic <= 0) then |
---|
705 | save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic) |
---|
706 | root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic) |
---|
707 | else |
---|
708 | save_var => Agrif_Curgrid % tabvars(indic) |
---|
709 | root_var => Agrif_Mygrid % tabvars(indic) |
---|
710 | endif |
---|
711 | ! |
---|
712 | if (.not.allocated(save_var%array4)) then |
---|
713 | allocate(save_var%array4(save_var%lb(1):save_var%ub(1),& |
---|
714 | save_var%lb(2):save_var%ub(2),& |
---|
715 | save_var%lb(3):save_var%ub(3),& |
---|
716 | save_var%lb(4):save_var%ub(4))) |
---|
717 | endif |
---|
718 | ! |
---|
719 | save_var % array4 = q |
---|
720 | root_var % restore = .true. |
---|
721 | !--------------------------------------------------------------------------------------------------- |
---|
722 | end subroutine Agrif_Save_ForRestore4D |
---|
723 | !=================================================================================================== |
---|
724 | ! |
---|
725 | end module Agrif_BcFunction |
---|