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 to initialize a fine grid from its parent grid, by using a space interpolation |
---|
25 | ! |
---|
26 | module Agrif_Interpolation |
---|
27 | ! |
---|
28 | use Agrif_InterpBasic |
---|
29 | use Agrif_Arrays |
---|
30 | use Agrif_Mask |
---|
31 | use Agrif_CurgridFunctions |
---|
32 | #if defined AGRIF_MPI |
---|
33 | use Agrif_Mpp |
---|
34 | #endif |
---|
35 | ! |
---|
36 | implicit none |
---|
37 | ! |
---|
38 | logical, private:: precomputedone(7) = .FALSE. |
---|
39 | ! |
---|
40 | private :: Agrif_Parentbounds |
---|
41 | private :: Agrif_Interp_1D_recursive, Agrif_Interp_2D_recursive, Agrif_Interp_3D_recursive |
---|
42 | private :: Agrif_Interp_4D_recursive, Agrif_Interp_5D_recursive, Agrif_Interp_6D_recursive |
---|
43 | private :: Agrif_InterpBase |
---|
44 | private :: Agrif_Find_list_interp, Agrif_AddTo_list_interp |
---|
45 | ! |
---|
46 | contains |
---|
47 | ! |
---|
48 | !=================================================================================================== |
---|
49 | ! subroutine Agrif_InterpVariable |
---|
50 | ! |
---|
51 | !> Sets some arguments of subroutine Agrif_InterpnD, n being the dimension of the grid variable |
---|
52 | !--------------------------------------------------------------------------------------------------- |
---|
53 | subroutine Agrif_InterpVariable ( parent, child, torestore, procname ) |
---|
54 | !--------------------------------------------------------------------------------------------------- |
---|
55 | type(Agrif_Variable), pointer :: parent !< Variable on the parent grid |
---|
56 | type(Agrif_Variable), pointer :: child !< Variable on the child grid |
---|
57 | logical, intent(in) :: torestore !< .false. indicates that the results of the |
---|
58 | !! interpolation are applied on the whole current grid |
---|
59 | procedure() :: procname !< Data recovery procedure |
---|
60 | !--------------------------------------------------------------------------------------------------- |
---|
61 | logical :: memberin |
---|
62 | integer :: nbdim ! Number of dimensions of the current grid |
---|
63 | integer, dimension(6) :: type_interp ! Type of interpolation (linear,spline,...) |
---|
64 | integer, dimension(6) :: nb_child |
---|
65 | integer, dimension(6) :: lb_child |
---|
66 | integer, dimension(6) :: ub_child |
---|
67 | integer, dimension(6) :: lb_parent |
---|
68 | real , dimension(6) :: s_child, s_parent |
---|
69 | real , dimension(6) :: ds_child, ds_parent |
---|
70 | integer, dimension(child % root_var % nbdim,2,2) :: childarray |
---|
71 | ! |
---|
72 | nbdim = child % root_var % nbdim |
---|
73 | type_interp = child % root_var % type_interp |
---|
74 | ! |
---|
75 | call PreProcessToInterpOrUpdate( parent, child, & |
---|
76 | nb_child, ub_child, & |
---|
77 | lb_child, lb_parent, & |
---|
78 | s_child, s_parent, & |
---|
79 | ds_child, ds_parent, nbdim, interp=.true.) |
---|
80 | ! |
---|
81 | ! Call to a procedure of interpolation against the number of dimensions of the grid variable |
---|
82 | ! |
---|
83 | call Agrif_InterpnD(type_interp, parent, child, & |
---|
84 | lb_child, ub_child, & |
---|
85 | lb_child, lb_parent, & |
---|
86 | s_child, s_parent, & |
---|
87 | ds_child, ds_parent, & |
---|
88 | child, torestore, nbdim, & |
---|
89 | childarray, memberin, & |
---|
90 | .false., procname, 0, 0) |
---|
91 | !--------------------------------------------------------------------------------------------------- |
---|
92 | end subroutine Agrif_InterpVariable |
---|
93 | !=================================================================================================== |
---|
94 | ! |
---|
95 | !=================================================================================================== |
---|
96 | ! subroutine Agrif_InterpnD |
---|
97 | ! |
---|
98 | !> Interpolates a nD grid variable from its parent grid, by using a space interpolation |
---|
99 | !--------------------------------------------------------------------------------------------------- |
---|
100 | subroutine Agrif_InterpnD ( type_interp, parent, child, pttab, petab, pttab_Child, pttab_Parent, & |
---|
101 | s_Child, s_Parent, ds_Child, ds_Parent, restore, torestore, & |
---|
102 | nbdim, childarray, memberin, in_bc, procname, nb, ndir ) |
---|
103 | !--------------------------------------------------------------------------------------------------- |
---|
104 | #if defined AGRIF_MPI |
---|
105 | include 'mpif.h' |
---|
106 | #endif |
---|
107 | ! |
---|
108 | INTEGER, DIMENSION(6), INTENT(in) :: type_interp !< Type of interpolation ! (linear,...) |
---|
109 | TYPE(Agrif_Variable), pointer :: parent !< Variable of the parent grid |
---|
110 | TYPE(Agrif_Variable), pointer :: child !< Variable of the child grid |
---|
111 | INTEGER, DIMENSION(nbdim), INTENT(in) :: pttab !< Index of the first point inside the domain |
---|
112 | INTEGER, DIMENSION(nbdim), INTENT(in) :: petab !< Index of the first point inside the domain |
---|
113 | INTEGER, DIMENSION(nbdim), INTENT(in) :: pttab_Child !< Index of the first point inside the domain |
---|
114 | !< for the child grid variable |
---|
115 | INTEGER, DIMENSION(nbdim), INTENT(in) :: pttab_Parent !< Index of the first point inside the domain |
---|
116 | !< for the parent grid variable |
---|
117 | REAL, DIMENSION(nbdim), INTENT(in) :: s_Child,s_Parent !< Positions of the parent and child grids |
---|
118 | REAL, DIMENSION(nbdim), INTENT(in) :: ds_Child,ds_Parent !< Space steps of the parent and child grids |
---|
119 | TYPE(Agrif_Variable), pointer :: restore !< Indicates points where interpolation |
---|
120 | LOGICAL, INTENT(in) :: torestore !< Indicates if the array restore is used |
---|
121 | INTEGER, INTENT(in) :: nbdim |
---|
122 | LOGICAL, INTENT(out) :: memberin |
---|
123 | LOGICAL, INTENT(in) :: in_bc !< .true. if called from Agrif_CorrectVariable \n |
---|
124 | !! .false. if called from Agrif_InterpVariable |
---|
125 | procedure() :: procname !< Data recovery procedure |
---|
126 | INTEGER, INTENT(in) :: nb, ndir |
---|
127 | ! |
---|
128 | INTEGER :: i,j,k,l,m,n |
---|
129 | integer :: i1,j1,k1 |
---|
130 | INTEGER, DIMENSION(nbdim) :: pttruetab,cetruetab |
---|
131 | INTEGER, DIMENSION(nbdim) :: indmin, indmax |
---|
132 | INTEGER, DIMENSION(nbdim) :: indminglob, indmaxglob |
---|
133 | #if defined AGRIF_MPI |
---|
134 | INTEGER, DIMENSION(nbdim) :: indminglob2,indmaxglob2 |
---|
135 | INTEGER, DIMENSION(nbdim) :: indminglob3,indmaxglob3 |
---|
136 | INTEGER, DIMENSION(:,:),ALLOCATABLE :: indminglob_chunks, indmaxglob_chunks |
---|
137 | INTEGER, DIMENSION(:,:),ALLOCATABLE :: indminglob2_chunks,indmaxglob2_chunks |
---|
138 | INTEGER, DIMENSION(:,:),ALLOCATABLE :: indminglob3_chunks,indmaxglob3_chunks |
---|
139 | #endif |
---|
140 | LOGICAL, DIMENSION(nbdim) :: noraftab |
---|
141 | REAL , DIMENSION(nbdim) :: s_Child_temp,s_Parent_temp |
---|
142 | INTEGER, DIMENSION(nbdim) :: lowerbound, upperbound, coords |
---|
143 | INTEGER, DIMENSION(nbdim,2,2), INTENT(OUT) :: childarray |
---|
144 | INTEGER, DIMENSION(nbdim,2,2) :: parentarray |
---|
145 | INTEGER, DIMENSION(nbdim,2,2) :: parentarray_decal |
---|
146 | LOGICAL :: member |
---|
147 | LOGICAL,DIMENSION(:),ALLOCATABLE :: member_chuncks |
---|
148 | INTEGER,DIMENSION(:,:),ALLOCATABLE :: decal_chunks |
---|
149 | LOGICAL :: find_list_interp |
---|
150 | ! |
---|
151 | #if defined AGRIF_MPI |
---|
152 | ! |
---|
153 | INTEGER, PARAMETER :: etiquette = 100 |
---|
154 | INTEGER :: code, local_proc |
---|
155 | INTEGER, DIMENSION(nbdim,4) :: tab3 |
---|
156 | INTEGER, DIMENSION(nbdim,4,0:Agrif_Nbprocs-1) :: tab4 |
---|
157 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8) :: tab4t |
---|
158 | INTEGER, DIMENSION(nbdim,2) :: tab5 |
---|
159 | INTEGER, DIMENSION(nbdim,2,0:Agrif_Nbprocs-1) :: tab6 |
---|
160 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,2) :: tab5t |
---|
161 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: memberinall |
---|
162 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: sendtoproc1, recvfromproc1 |
---|
163 | LOGICAL, DIMENSION(1) :: memberin1 |
---|
164 | LOGICAL :: memberout |
---|
165 | ! |
---|
166 | #endif |
---|
167 | ! CHUNK (periodicity) |
---|
168 | INTEGER :: nb_chunks |
---|
169 | INTEGER :: agrif_external_switch_index |
---|
170 | INTEGER, DIMENSION(2) :: test_orientation |
---|
171 | !INTEGER, DIMENSION(2,nbdim,2,2) :: parentarray_chunk |
---|
172 | INTEGER, DIMENSION(:,:,:,:), ALLOCATABLE :: parentarray_chunk |
---|
173 | INTEGER, DIMENSION(:,:,:,:), ALLOCATABLE :: parentarray_chunk_decal |
---|
174 | INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: bounds_chunks |
---|
175 | logical,dimension(:),allocatable :: correction_required |
---|
176 | ! |
---|
177 | type(Agrif_Variable), pointer, save :: tempC => NULL() ! Temporary child grid variable |
---|
178 | type(Agrif_Variable), pointer, save :: tempP => NULL() ! Temporary parent grid variable |
---|
179 | type(Agrif_Variable), pointer, save :: tempPextend => NULL() ! Temporary parent grid variable |
---|
180 | type(Agrif_Variable), pointer, save :: parentvalues => NULL() |
---|
181 | ! |
---|
182 | coords = child % root_var % coords |
---|
183 | ! |
---|
184 | ! Boundaries of the current grid where interpolation is done |
---|
185 | find_list_interp = & |
---|
186 | Agrif_Find_list_interp( & |
---|
187 | child % list_interp, & |
---|
188 | pttab, petab, pttab_Child, pttab_Parent, nbdim, & |
---|
189 | indmin, indmax, indminglob, indmaxglob, & |
---|
190 | pttruetab, cetruetab, memberin, & |
---|
191 | parentarray_chunk,parentarray_chunk_decal,decal_chunks, & |
---|
192 | correction_required,member_chuncks,nb_chunks & |
---|
193 | #if defined AGRIF_MPI |
---|
194 | ,indminglob2, indmaxglob2, parentarray, & |
---|
195 | member, tab4t,memberinall, sendtoproc1, recvfromproc1 & |
---|
196 | #endif |
---|
197 | ) |
---|
198 | ! |
---|
199 | if (.not.find_list_interp) then |
---|
200 | ! |
---|
201 | call Agrif_get_var_bounds_array(child, lowerbound, upperbound, nbdim, parent) |
---|
202 | call Agrif_Childbounds(nbdim, lowerbound, upperbound, & |
---|
203 | pttab, petab, Agrif_Procrank, coords, & |
---|
204 | pttruetab, cetruetab, memberin) |
---|
205 | |
---|
206 | if (agrif_debug_interp) then |
---|
207 | print *,'************CHILDBOUNDS*********************************' |
---|
208 | #ifdef AGRIF_MPI |
---|
209 | print *,'Processeur ',Agrif_Procrank |
---|
210 | #endif |
---|
211 | print *,'memberin ',memberin |
---|
212 | do i = 1 , nbdim |
---|
213 | print *,'Direction ',i,' indices debut: ',pttab(i),pttruetab(i) |
---|
214 | print *,'Direction ',i,' indices fin : ',petab(i),cetruetab(i) |
---|
215 | enddo |
---|
216 | print *,'*********************************************' |
---|
217 | endif |
---|
218 | |
---|
219 | call Agrif_Parentbounds(type_interp,nbdim,indminglob,indmaxglob, & |
---|
220 | s_Parent_temp,s_Child_temp, & |
---|
221 | s_Child,ds_Child, & |
---|
222 | s_Parent,ds_Parent, & |
---|
223 | pttab,petab, & |
---|
224 | pttab_Child,pttab_Parent, & |
---|
225 | child%root_var % posvar, coords) |
---|
226 | |
---|
227 | if (agrif_debug_interp) then |
---|
228 | print *,'************PARENTBOUNDS*********************************' |
---|
229 | #ifdef AGRIF_MPI |
---|
230 | print *,'Processeur ',Agrif_Procrank |
---|
231 | #endif |
---|
232 | do i = 1 , nbdim |
---|
233 | print *,'Direction ',i,' indices debut: ',pttab(i),indminglob(i) |
---|
234 | print *,'Direction ',i,' indices fin : ',petab(i),indmaxglob(i) |
---|
235 | enddo |
---|
236 | |
---|
237 | do i = 1 , nbdim |
---|
238 | print *,'Direction ',i,' s_parent_temp: ',s_parent_temp(i) |
---|
239 | print *,'Direction ',i,' s_Child_temp : ',s_Child_temp(i) |
---|
240 | enddo |
---|
241 | print *,'*********************************************' |
---|
242 | endif |
---|
243 | |
---|
244 | #if defined AGRIF_MPI |
---|
245 | if (memberin) then |
---|
246 | call Agrif_Parentbounds(type_interp,nbdim,indmin,indmax, & |
---|
247 | s_Parent_temp,s_Child_temp, & |
---|
248 | s_Child,ds_Child, & |
---|
249 | s_Parent,ds_Parent, & |
---|
250 | pttruetab,cetruetab, & |
---|
251 | pttab_Child,pttab_Parent, & |
---|
252 | child%root_var % posvar, coords) |
---|
253 | |
---|
254 | endif |
---|
255 | if (agrif_debug_interp) then |
---|
256 | print *,'************PARENTBOUNDSMPI*********************************' |
---|
257 | #ifdef AGRIF_MPI |
---|
258 | print *,'Processeur ',Agrif_Procrank |
---|
259 | #endif |
---|
260 | do i = 1 , nbdim |
---|
261 | print *,'Direction ',i,' indices debut: ',pttruetab(i),indmin(i) |
---|
262 | print *,'Direction ',i,' indices fin : ',cetruetab(i),indmax(i) |
---|
263 | enddo |
---|
264 | |
---|
265 | do i = 1 , nbdim |
---|
266 | print *,'Direction ',i,' s_parent_temp: ',s_parent_temp(i) |
---|
267 | print *,'Direction ',i,' s_Child_temp : ',s_Child_temp(i) |
---|
268 | enddo |
---|
269 | print *,'*********************************************' |
---|
270 | endif |
---|
271 | |
---|
272 | local_proc = Agrif_Procrank |
---|
273 | call Agrif_get_var_bounds_array(parent,lowerbound,upperbound,nbdim) |
---|
274 | call Agrif_ChildGrid_to_ParentGrid() |
---|
275 | |
---|
276 | parentarray(:,1,1) = indminglob |
---|
277 | parentarray(:,2,1) = indmaxglob |
---|
278 | parentarray(:,1,2) = indminglob |
---|
279 | parentarray(:,2,2) = indmaxglob |
---|
280 | |
---|
281 | if (associated(agrif_external_mapping)) then |
---|
282 | |
---|
283 | call agrif_external_mapping(nbdim,child%root_var % posvar(1),child%root_var % posvar(2), & |
---|
284 | parentarray,parentarray_chunk,correction_required,nb_chunks) |
---|
285 | allocate(decal_chunks(nb_chunks,nbdim)) |
---|
286 | do i=1,nb_chunks |
---|
287 | decal_chunks(i,:)=parentarray_chunk(i,:,1,1)-parentarray_chunk(i,:,1,2) |
---|
288 | enddo |
---|
289 | else |
---|
290 | nb_chunks=1 |
---|
291 | allocate(correction_required(nb_chunks)) |
---|
292 | correction_required=.FALSE. |
---|
293 | allocate(parentarray_chunk(nb_chunks,nbdim,2,2)) |
---|
294 | parentarray_chunk(1,:,:,:)=parentarray |
---|
295 | allocate(decal_chunks(nb_chunks,nbdim)) |
---|
296 | decal_chunks=0 |
---|
297 | endif |
---|
298 | if (agrif_debug_interp) then |
---|
299 | print *,'AVANT PARENTCHILDBOUNDS' |
---|
300 | print *,'nombre de chunks ',nb_chunks |
---|
301 | do i=1,nb_chunks |
---|
302 | print *,'CHUNK Number : ',i |
---|
303 | do j=1,nbdim |
---|
304 | print *,'Direction ',j |
---|
305 | print *,'MIN MAX (2) = ',parentarray_chunk(i,j,1,2),parentarray_chunk(i,j,2,2) |
---|
306 | print *,'MIN MAX (1) = ',parentarray_chunk(i,j,1,1),parentarray_chunk(i,j,2,1) |
---|
307 | enddo |
---|
308 | enddo |
---|
309 | print *,'APRES PARENTCHILDBOUNDS' |
---|
310 | endif |
---|
311 | |
---|
312 | allocate(indminglob_chunks(nb_chunks,nbdim)) |
---|
313 | allocate(indmaxglob_chunks(nb_chunks,nbdim)) |
---|
314 | allocate(indminglob2_chunks(nb_chunks,nbdim)) |
---|
315 | allocate(indmaxglob2_chunks(nb_chunks,nbdim)) |
---|
316 | allocate(indminglob3_chunks(nb_chunks,nbdim)) |
---|
317 | allocate(indmaxglob3_chunks(nb_chunks,nbdim)) |
---|
318 | allocate(member_chuncks(nb_chunks)) |
---|
319 | |
---|
320 | do i=1,nb_chunks |
---|
321 | indminglob_chunks(i,:) = parentarray_chunk(i,:,1,2) |
---|
322 | indmaxglob_chunks(i,:) = parentarray_chunk(i,:,2,2) |
---|
323 | enddo |
---|
324 | |
---|
325 | do i=1,nb_chunks |
---|
326 | |
---|
327 | call Agrif_Childbounds(nbdim,lowerbound,upperbound, & |
---|
328 | indminglob_chunks(i,:),indmaxglob_chunks(i,:), local_proc, coords, & |
---|
329 | indminglob2_chunks(i,:),indmaxglob2_chunks(i,:),member_chuncks(i), & |
---|
330 | indminglob3_chunks(i,:),indmaxglob3_chunks(i,:)) |
---|
331 | enddo |
---|
332 | ! |
---|
333 | ! call Agrif_Childbounds(nbdim,lowerbound,upperbound, & |
---|
334 | ! indminglob,indmaxglob, local_proc, coords, & |
---|
335 | ! indminglob2,indmaxglob2,member, & |
---|
336 | ! indminglob3,indmaxglob3,check_perio=.TRUE.) |
---|
337 | |
---|
338 | if (agrif_debug_interp) then |
---|
339 | print *,'************CHILDBOUNDSPARENTMPI*********************************' |
---|
340 | #ifdef AGRIF_MPI |
---|
341 | print *,'Processeur ',Agrif_Procrank |
---|
342 | #endif |
---|
343 | do j=1,nb_chunks |
---|
344 | print *,'Chunk number ',j |
---|
345 | |
---|
346 | do i = 1 , nbdim |
---|
347 | print *,'Direction ',i,' indices debut: ',indminglob_chunks(j,i),indminglob2_chunks(j,i),indminglob3_chunks(j,i) |
---|
348 | print *,'Direction ',i,' indices fin : ',indmaxglob_chunks(j,i),indmaxglob2_chunks(j,i),indmaxglob3_chunks(j,i) |
---|
349 | enddo |
---|
350 | enddo |
---|
351 | print *,'*********************************************' |
---|
352 | endif |
---|
353 | ! |
---|
354 | ! if (member) then |
---|
355 | ! call Agrif_GlobalToLocalBounds(parentarray, & |
---|
356 | ! lowerbound, upperbound, & |
---|
357 | ! indminglob2, indmaxglob2, coords,& |
---|
358 | ! nbdim, local_proc, member,check_perio=.TRUE.) |
---|
359 | ! if (agrif_debug_interp) then |
---|
360 | ! do i=1,nbdim |
---|
361 | ! print *,'parentarray = ',i,parentarray(i,1,1),parentarray(i,2,1), & |
---|
362 | ! parentarray(i,1,2),parentarray(i,2,2) |
---|
363 | ! enddo |
---|
364 | ! endif |
---|
365 | ! endif |
---|
366 | |
---|
367 | allocate(parentarray_chunk_decal(nb_chunks,nbdim,2,2)) |
---|
368 | do j=1,nb_chunks |
---|
369 | if (agrif_debug_interp) print *,'CHUNK = ',j |
---|
370 | if (member_chuncks(j)) then |
---|
371 | ! call Agrif_GlobalToLocalBounds(parentarray_chunk(j,:,:,:), & |
---|
372 | ! lowerbound, upperbound, & |
---|
373 | ! indminglob2_chunks(j,:), indmaxglob2_chunks(j,:), coords, & |
---|
374 | ! nbdim, local_proc, member_chuncks(j),check_perio=.TRUE.) |
---|
375 | |
---|
376 | call Agrif_GlobalToLocalBounds(parentarray_chunk(j,:,:,:), & |
---|
377 | lowerbound, upperbound, & |
---|
378 | indminglob2_chunks(j,:), indmaxglob2_chunks(j,:), coords, & |
---|
379 | nbdim, local_proc, member_chuncks(j)) |
---|
380 | |
---|
381 | if (correction_required(j)) then |
---|
382 | do i=1,2 |
---|
383 | test_orientation(1)=agrif_external_switch_index(child%root_var % posvar(1),child%root_var % posvar(2), & |
---|
384 | parentarray_chunk(j,i,1,1),i) |
---|
385 | test_orientation(2)=agrif_external_switch_index(child%root_var % posvar(1),child%root_var % posvar(2), & |
---|
386 | parentarray_chunk(j,i,2,1),i) |
---|
387 | parentarray_chunk_decal(j,i,1,1)=minval(test_orientation) |
---|
388 | parentarray_chunk_decal(j,i,2,1)=maxval(test_orientation) |
---|
389 | enddo |
---|
390 | do i=3,nbdim |
---|
391 | parentarray_chunk_decal(j,i,:,1) = parentarray_chunk(j,i,:,1)+decal_chunks(j,i) |
---|
392 | enddo |
---|
393 | else |
---|
394 | do i=1,nbdim |
---|
395 | parentarray_chunk_decal(j,i,:,1) = parentarray_chunk(j,i,:,1)+decal_chunks(j,i) |
---|
396 | enddo |
---|
397 | endif |
---|
398 | |
---|
399 | if (agrif_debug_interp) then |
---|
400 | do i=1,nbdim |
---|
401 | print *,'parentarray = ',i,parentarray_chunk(j,i,1,1),parentarray_chunk(j,i,2,1), & |
---|
402 | parentarray_chunk(j,i,1,2),parentarray_chunk(j,i,2,2) |
---|
403 | print *,'parentarraydecal = ',i,parentarray_chunk_decal(j,i,1,1),parentarray_chunk_decal(j,i,2,1) |
---|
404 | enddo |
---|
405 | endif |
---|
406 | endif |
---|
407 | enddo |
---|
408 | |
---|
409 | parentarray(:,1,:)=Huge(1) |
---|
410 | parentarray(:,2,:)=-Huge(1) |
---|
411 | indminglob2=Huge(1) |
---|
412 | indmaxglob2=-Huge(1) |
---|
413 | indminglob3=Huge(1) |
---|
414 | indmaxglob3=-Huge(1) |
---|
415 | member = .FALSE. |
---|
416 | do j=1,nb_chunks |
---|
417 | if (member_chuncks(j)) then |
---|
418 | do i=1,nbdim |
---|
419 | parentarray(i,1,1) = min(parentarray(i,1,1),parentarray_chunk_decal(j,i,1,1)) |
---|
420 | parentarray(i,1,2) = min(parentarray(i,1,2),parentarray_chunk(j,i,1,2)) |
---|
421 | parentarray(i,2,1) = max(parentarray(i,2,1),parentarray_chunk_decal(j,i,2,1)) |
---|
422 | parentarray(i,2,2) = max(parentarray(i,2,2),parentarray_chunk(j,i,2,2)) |
---|
423 | enddo |
---|
424 | if (correction_required(j)) then |
---|
425 | if (agrif_debug_interp) then |
---|
426 | do i=1,nbdim |
---|
427 | print *,'direction ',i |
---|
428 | print *,'glob2_chuk = ',indminglob2_chunks(j,i),indmaxglob2_chunks(j,i) |
---|
429 | print *,'glob3_chuk = ',indminglob3_chunks(j,i),indmaxglob3_chunks(j,i) |
---|
430 | enddo |
---|
431 | endif |
---|
432 | do i=1,2 |
---|
433 | test_orientation(1)=agrif_external_switch_index(child%root_var % posvar(1),child%root_var % posvar(2), & |
---|
434 | indminglob2_chunks(j,i),i) |
---|
435 | test_orientation(2)=agrif_external_switch_index(child%root_var % posvar(1),child%root_var % posvar(2), & |
---|
436 | indmaxglob2_chunks(j,i),i) |
---|
437 | indminglob2(i)=min(indminglob2(i),minval(test_orientation)) |
---|
438 | indmaxglob2(i)=max(indmaxglob2(i),maxval(test_orientation)) |
---|
439 | enddo |
---|
440 | |
---|
441 | do i=1,2 |
---|
442 | test_orientation(1)=agrif_external_switch_index(child%root_var % posvar(1),child%root_var % posvar(2), & |
---|
443 | indminglob3_chunks(j,i),i) |
---|
444 | test_orientation(2)=agrif_external_switch_index(child%root_var % posvar(1),child%root_var % posvar(2), & |
---|
445 | indmaxglob3_chunks(j,i),i) |
---|
446 | indminglob3(i)=min(indminglob3(i),minval(test_orientation)) |
---|
447 | indmaxglob3(i)=max(indmaxglob3(i),maxval(test_orientation)) |
---|
448 | enddo |
---|
449 | |
---|
450 | do i=3,nbdim |
---|
451 | indminglob2(i)=min(indminglob2(i),indminglob2_chunks(j,i)+decal_chunks(j,i)) |
---|
452 | indmaxglob2(i)=max(indmaxglob2(i),indmaxglob2_chunks(j,i)+decal_chunks(j,i)) |
---|
453 | indminglob3(i)=min(indminglob3(i),indminglob3_chunks(j,i)+decal_chunks(j,i)) |
---|
454 | indmaxglob3(i)=max(indmaxglob3(i),indmaxglob3_chunks(j,i)+decal_chunks(j,i)) |
---|
455 | enddo |
---|
456 | else |
---|
457 | do i=1,nbdim |
---|
458 | indminglob2(i)=min(indminglob2(i),indminglob2_chunks(j,i)+decal_chunks(j,i)) |
---|
459 | indmaxglob2(i)=max(indmaxglob2(i),indmaxglob2_chunks(j,i)+decal_chunks(j,i)) |
---|
460 | indminglob3(i)=min(indminglob3(i),indminglob3_chunks(j,i)+decal_chunks(j,i)) |
---|
461 | indmaxglob3(i)=max(indmaxglob3(i),indmaxglob3_chunks(j,i)+decal_chunks(j,i)) |
---|
462 | enddo |
---|
463 | endif |
---|
464 | |
---|
465 | member = .TRUE. |
---|
466 | endif |
---|
467 | enddo |
---|
468 | |
---|
469 | call Agrif_ParentGrid_to_ChildGrid() |
---|
470 | |
---|
471 | if (agrif_debug_interp) then |
---|
472 | print *,'************ FINAL PARENTARRAY *****************' |
---|
473 | #ifdef AGRIF_MPI |
---|
474 | print *,'Processeur ',Agrif_Procrank,' MEMBER = ',member |
---|
475 | do i=1,nbdim |
---|
476 | print *,'Direction ',i,' indices debut = ',parentarray(i,1,1),parentarray(i,1,2) |
---|
477 | print *,'Direction ',i,' indices fin = ',parentarray(i,2,1),parentarray(i,2,2) |
---|
478 | enddo |
---|
479 | #endif |
---|
480 | endif |
---|
481 | |
---|
482 | if (agrif_debug_interp) then |
---|
483 | print *,'************ FINAL INDMINGLOB *****************' |
---|
484 | #ifdef AGRIF_MPI |
---|
485 | print *,'Processeur ',Agrif_Procrank,' MEMBER = ',member |
---|
486 | do i=1,nbdim |
---|
487 | print *,'Direction ',i,' indices debut = ',indminglob2(i),indminglob3(i) |
---|
488 | print *,'Direction ',i,' indices fin = ',indmaxglob2(i),indmaxglob3(i) |
---|
489 | enddo |
---|
490 | #endif |
---|
491 | endif |
---|
492 | |
---|
493 | #else |
---|
494 | parentarray(:,1,1) = indminglob |
---|
495 | parentarray(:,2,1) = indmaxglob |
---|
496 | parentarray(:,1,2) = indminglob |
---|
497 | parentarray(:,2,2) = indmaxglob |
---|
498 | |
---|
499 | |
---|
500 | if (associated(agrif_external_mapping)) then |
---|
501 | call Agrif_ChildGrid_to_ParentGrid() |
---|
502 | call agrif_external_mapping(nbdim,child%root_var % posvar(1),child%root_var % posvar(2), & |
---|
503 | parentarray,parentarray_chunk,correction_required,nb_chunks) |
---|
504 | call Agrif_ParentGrid_to_ChildGrid() |
---|
505 | allocate(decal_chunks(nb_chunks,nbdim)) |
---|
506 | do i=1,nb_chunks |
---|
507 | decal_chunks(i,:)=parentarray_chunk(i,:,1,1)-parentarray_chunk(i,:,1,2) |
---|
508 | enddo |
---|
509 | else |
---|
510 | nb_chunks=1 |
---|
511 | allocate(correction_required(nb_chunks)) |
---|
512 | correction_required=.FALSE. |
---|
513 | allocate(parentarray_chunk(nb_chunks,nbdim,2,2)) |
---|
514 | parentarray_chunk(1,:,:,:)=parentarray |
---|
515 | endif |
---|
516 | if (agrif_debug_interp) then |
---|
517 | print *,'AVANT PARENTCHILDBOUNDS' |
---|
518 | print *,'nombre de chunks ',nb_chunks |
---|
519 | do i=1,nb_chunks |
---|
520 | print *,'CHUNK Number : ',i |
---|
521 | do j=1,nbdim |
---|
522 | print *,'Direction ',j |
---|
523 | print *,'MIN MAX (2) = ',parentarray_chunk(i,j,1,2),parentarray_chunk(i,j,2,2) |
---|
524 | print *,'MIN MAX (1) = ',parentarray_chunk(i,j,1,1),parentarray_chunk(i,j,2,1) |
---|
525 | enddo |
---|
526 | enddo |
---|
527 | print *,'APRES PARENTCHILDBOUNDS' |
---|
528 | endif |
---|
529 | allocate(member_chuncks(nb_chunks)) |
---|
530 | allocate(parentarray_chunk_decal(nb_chunks,nbdim,2,2)) |
---|
531 | member_chuncks = .TRUE. |
---|
532 | member = .TRUE. |
---|
533 | do j=1,nb_chunks |
---|
534 | if (agrif_debug_interp) print *,'CHUNK = ',j |
---|
535 | if (member_chuncks(j)) then |
---|
536 | do i=1,nbdim |
---|
537 | parentarray_chunk_decal(j,i,:,1) = parentarray_chunk(j,i,:,1) !+decal_chunks(j,i) |
---|
538 | if (agrif_debug_interp) then |
---|
539 | print *,'ENCORE = ',parentarray_chunk(j,i,:,1),parentarray_chunk_decal(j,i,:,1) |
---|
540 | endif |
---|
541 | enddo |
---|
542 | if (agrif_debug_interp) then |
---|
543 | do i=1,nbdim |
---|
544 | print *,'parentarray = ',i,parentarray_chunk(j,i,1,1),parentarray_chunk(j,i,2,1), & |
---|
545 | parentarray_chunk(j,i,1,2),parentarray_chunk(j,i,2,2) |
---|
546 | print *,'parentarraydecal = ',i,parentarray_chunk_decal(j,i,1,1),parentarray_chunk_decal(j,i,2,1) |
---|
547 | enddo |
---|
548 | endif |
---|
549 | endif |
---|
550 | enddo |
---|
551 | |
---|
552 | |
---|
553 | indmin = indminglob |
---|
554 | indmax = indmaxglob |
---|
555 | member = .TRUE. |
---|
556 | #endif |
---|
557 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
558 | ! Correct for non refined directions |
---|
559 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
560 | do i=1,nbdim |
---|
561 | if (coords(i) == 0) then |
---|
562 | indmin(i) = indminglob(i) |
---|
563 | indmax(i) = indmaxglob(i) |
---|
564 | pttruetab(i) = indminglob(i) |
---|
565 | cetruetab(i) = indmaxglob(i) |
---|
566 | endif |
---|
567 | enddo |
---|
568 | |
---|
569 | else |
---|
570 | |
---|
571 | #if defined AGRIF_MPI |
---|
572 | s_Parent_temp = s_Parent + (indmin - pttab_Parent) * ds_Parent |
---|
573 | s_Child_temp = s_Child + (pttruetab - pttab_Child) * ds_Child |
---|
574 | #else |
---|
575 | parentarray(:,1,1) = indminglob |
---|
576 | parentarray(:,2,1) = indmaxglob |
---|
577 | parentarray(:,1,2) = indminglob |
---|
578 | parentarray(:,2,2) = indmaxglob |
---|
579 | indmin = indminglob |
---|
580 | indmax = indmaxglob |
---|
581 | member = .TRUE. |
---|
582 | s_Parent_temp = s_Parent + (indminglob - pttab_Parent) * ds_Parent |
---|
583 | s_Child_temp = s_Child + (pttab - pttab_Child) * ds_Child |
---|
584 | #endif |
---|
585 | |
---|
586 | endif |
---|
587 | |
---|
588 | if (agrif_debug_interp) then |
---|
589 | print *,'************SPARENTCHILD*********************************' |
---|
590 | #ifdef AGRIF_MPI |
---|
591 | print *,'Processeur ',Agrif_Procrank |
---|
592 | #endif |
---|
593 | do i = 1 , nbdim |
---|
594 | print *,'Direction ',i,' s_Parent_temp : ',s_Parent_temp(i),indmin(i) |
---|
595 | print *,'Direction ',i,' s_Child_temp : ',s_Child_temp(i),pttruetab(i) |
---|
596 | enddo |
---|
597 | print *,'*********************************************' |
---|
598 | endif |
---|
599 | |
---|
600 | call Agrif_get_var_bounds_array(parent,lowerbound,upperbound,nbdim) |
---|
601 | ! |
---|
602 | if (member) then |
---|
603 | if (.not.associated(tempP)) allocate(tempP) |
---|
604 | ! |
---|
605 | call Agrif_array_allocate(tempP,parentarray(:,1,1),parentarray(:,2,1),nbdim) |
---|
606 | call Agrif_var_set_array_tozero(tempP,nbdim) |
---|
607 | endif |
---|
608 | Agrif_CurChildgrid=>Agrif_Curgrid |
---|
609 | call Agrif_ChildGrid_to_ParentGrid() |
---|
610 | do i=1,nb_chunks |
---|
611 | if (agrif_debug_interp) then |
---|
612 | print *,'PROCNAME POUR CHUCNK ',i |
---|
613 | endif |
---|
614 | |
---|
615 | if (member_chuncks(i)) then |
---|
616 | select case (nbdim) |
---|
617 | case(1) |
---|
618 | ! call procname(tempP%array1, & |
---|
619 | ! parentarray(1,1,2),parentarray(1,2,2),.TRUE.,nb,ndir) |
---|
620 | |
---|
621 | call procname(tempP%array1(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1)), & |
---|
622 | parentarray_chunk(i,1,1,2),parentarray_chunk(i,1,2,2),.TRUE.,nb,ndir) |
---|
623 | |
---|
624 | case(2) |
---|
625 | ! call procname(tempP%array2, & |
---|
626 | ! parentarray(1,1,2),parentarray(1,2,2), & |
---|
627 | ! parentarray(2,1,2),parentarray(2,2,2),.TRUE.,nb,ndir) |
---|
628 | |
---|
629 | call procname(tempP%array2(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
630 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1)), & |
---|
631 | parentarray_chunk(i,1,1,2),parentarray_chunk(i,1,2,2), & |
---|
632 | parentarray_chunk(i,2,1,2),parentarray_chunk(i,2,2,2),.TRUE.,nb,ndir) |
---|
633 | if (agrif_debug_interp) print *,'SORTIE DE PROCNAME' |
---|
634 | if (correction_required(i)) then |
---|
635 | call correct_field(tempP%array2(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
636 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1)), & |
---|
637 | parentarray_chunk_decal(i,1,1,1),parentarray_chunk_decal(i,1,2,1), & |
---|
638 | parentarray_chunk_decal(i,2,1,1),parentarray_chunk_decal(i,2,2,1)) |
---|
639 | endif |
---|
640 | |
---|
641 | case(3) |
---|
642 | call procname(tempP%array3(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
643 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1), & |
---|
644 | parentarray_chunk_decal(i,3,1,1):parentarray_chunk_decal(i,3,2,1)), & |
---|
645 | parentarray_chunk(i,1,1,2),parentarray_chunk(i,1,2,2), & |
---|
646 | parentarray_chunk(i,2,1,2),parentarray_chunk(i,2,2,2), & |
---|
647 | parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2),.TRUE.,nb,ndir) |
---|
648 | |
---|
649 | if (agrif_debug_interp) then |
---|
650 | print *,'CHUNK = ',i |
---|
651 | print *,'NBNDIR = ',nb,ndir,correction_required(i) |
---|
652 | print *,'TEMPARRAY3 INDEX LOCAUX PUIS GLOBAUX' |
---|
653 | print *,parentarray_chunk_decal(i,1,1,1),parentarray_chunk_decal(i,1,2,1) |
---|
654 | print *,parentarray_chunk_decal(i,2,1,1),parentarray_chunk_decal(i,2,2,1) |
---|
655 | print *,parentarray_chunk_decal(i,3,1,1),parentarray_chunk_decal(i,3,2,1) |
---|
656 | print *,parentarray_chunk(i,1,1,2),parentarray_chunk(i,1,2,2) |
---|
657 | print *,parentarray_chunk(i,2,1,2),parentarray_chunk(i,2,2,2) |
---|
658 | print *,parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2) |
---|
659 | do j1=parentarray_chunk_decal(i,2,1,1),parentarray_chunk_decal(i,2,2,1) |
---|
660 | do i1=parentarray_chunk_decal(i,1,1,1),parentarray_chunk_decal(i,1,2,1) |
---|
661 | print *,'valprocname = ',i1,j1,tempP%array3(i1,j1,1) |
---|
662 | enddo |
---|
663 | enddo |
---|
664 | endif |
---|
665 | if (correction_required(i)) then |
---|
666 | do k=parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2) |
---|
667 | call correct_field(tempP%array3(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
668 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1),k), & |
---|
669 | parentarray_chunk_decal(i,1,1,1),parentarray_chunk_decal(i,1,2,1), & |
---|
670 | parentarray_chunk_decal(i,2,1,1),parentarray_chunk_decal(i,2,2,1)) |
---|
671 | enddo |
---|
672 | if (agrif_debug_interp) then |
---|
673 | do j1=parentarray_chunk_decal(i,2,1,1),parentarray_chunk_decal(i,2,2,1) |
---|
674 | do i1=parentarray_chunk_decal(i,1,1,1),parentarray_chunk_decal(i,1,2,1) |
---|
675 | print *,'valprocname apres correction = ',i1,j1,tempP%array3(i1,j1,1) |
---|
676 | enddo |
---|
677 | enddo |
---|
678 | endif |
---|
679 | endif |
---|
680 | |
---|
681 | ! call procname(tempP%array3, & |
---|
682 | ! parentarray_chunk(i,1,1,2),parentarray_chunk(i,1,2,2), & |
---|
683 | ! parentarray_chunk(i,2,1,2),parentarray_chunk(i,2,2,2), & |
---|
684 | ! parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2),.TRUE.,nb,ndir) |
---|
685 | case(4) |
---|
686 | |
---|
687 | call procname(tempP%array4(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
688 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1), & |
---|
689 | parentarray_chunk_decal(i,3,1,1):parentarray_chunk_decal(i,3,2,1), & |
---|
690 | parentarray_chunk_decal(i,4,1,1):parentarray_chunk_decal(i,4,2,1)), & |
---|
691 | parentarray_chunk(i,1,1,2),parentarray_chunk(i,1,2,2), & |
---|
692 | parentarray_chunk(i,2,1,2),parentarray_chunk(i,2,2,2), & |
---|
693 | parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2), & |
---|
694 | parentarray_chunk(i,4,1,2),parentarray_chunk(i,4,2,2),.TRUE.,nb,ndir) |
---|
695 | |
---|
696 | if (correction_required(i)) then |
---|
697 | do l=parentarray_chunk(i,4,1,2),parentarray_chunk(i,4,2,2) |
---|
698 | do k=parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2) |
---|
699 | call correct_field(tempP%array4(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
700 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1),k,l), & |
---|
701 | parentarray_chunk_decal(i,1,1,1),parentarray_chunk_decal(i,1,2,1), & |
---|
702 | parentarray_chunk_decal(i,2,1,1),parentarray_chunk_decal(i,2,2,1)) |
---|
703 | enddo |
---|
704 | enddo |
---|
705 | endif |
---|
706 | |
---|
707 | ! call procname(tempP%array4, & |
---|
708 | ! parentarray(1,1,2),parentarray(1,2,2), & |
---|
709 | ! parentarray(2,1,2),parentarray(2,2,2), & |
---|
710 | ! parentarray(3,1,2),parentarray(3,2,2), & |
---|
711 | ! parentarray(4,1,2),parentarray(4,2,2),.TRUE.,nb,ndir) |
---|
712 | case(5) |
---|
713 | |
---|
714 | call procname(tempP%array5(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
715 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1), & |
---|
716 | parentarray_chunk_decal(i,3,1,1):parentarray_chunk_decal(i,3,2,1), & |
---|
717 | parentarray_chunk_decal(i,4,1,1):parentarray_chunk_decal(i,4,2,1), & |
---|
718 | parentarray_chunk_decal(i,5,1,1):parentarray_chunk_decal(i,5,2,1)), & |
---|
719 | parentarray_chunk(i,1,1,2),parentarray_chunk(i,1,2,2), & |
---|
720 | parentarray_chunk(i,2,1,2),parentarray_chunk(i,2,2,2), & |
---|
721 | parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2), & |
---|
722 | parentarray_chunk(i,4,1,2),parentarray_chunk(i,4,2,2), & |
---|
723 | parentarray_chunk(i,5,1,2),parentarray_chunk(i,5,2,2),.TRUE.,nb,ndir) |
---|
724 | |
---|
725 | if (correction_required(i)) then |
---|
726 | do m=parentarray_chunk(i,5,1,2),parentarray_chunk(i,5,2,2) |
---|
727 | do l=parentarray_chunk(i,4,1,2),parentarray_chunk(i,4,2,2) |
---|
728 | do k=parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2) |
---|
729 | call correct_field(tempP%array5(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
730 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1),k,l,m), & |
---|
731 | parentarray_chunk_decal(i,1,1,1),parentarray_chunk_decal(i,1,2,1), & |
---|
732 | parentarray_chunk_decal(i,2,1,1),parentarray_chunk_decal(i,2,2,1)) |
---|
733 | enddo |
---|
734 | enddo |
---|
735 | enddo |
---|
736 | endif |
---|
737 | |
---|
738 | ! call procname(tempP%array5, & |
---|
739 | ! parentarray(1,1,2),parentarray(1,2,2), & |
---|
740 | ! parentarray(2,1,2),parentarray(2,2,2), & |
---|
741 | ! parentarray(3,1,2),parentarray(3,2,2), & |
---|
742 | ! parentarray(4,1,2),parentarray(4,2,2), & |
---|
743 | ! parentarray(5,1,2),parentarray(5,2,2),.TRUE.,nb,ndir) |
---|
744 | case(6) |
---|
745 | |
---|
746 | call procname(tempP%array6(parentarray_chunk_decal(i,1,1,1):parentarray_chunk_decal(i,1,2,1), & |
---|
747 | parentarray_chunk_decal(i,2,1,1):parentarray_chunk_decal(i,2,2,1), & |
---|
748 | parentarray_chunk_decal(i,3,1,1):parentarray_chunk_decal(i,3,2,1), & |
---|
749 | parentarray_chunk_decal(i,4,1,1):parentarray_chunk_decal(i,4,2,1), & |
---|
750 | parentarray_chunk_decal(i,5,1,1):parentarray_chunk_decal(i,5,2,1), & |
---|
751 | parentarray_chunk_decal(i,6,1,1):parentarray_chunk_decal(i,6,2,1)), & |
---|
752 | parentarray_chunk(i,1,1,2),parentarray_chunk(i,1,2,2), & |
---|
753 | parentarray_chunk(i,2,1,2),parentarray_chunk(i,2,2,2), & |
---|
754 | parentarray_chunk(i,3,1,2),parentarray_chunk(i,3,2,2), & |
---|
755 | parentarray_chunk(i,4,1,2),parentarray_chunk(i,4,2,2), & |
---|
756 | parentarray_chunk(i,5,1,2),parentarray_chunk(i,5,2,2), & |
---|
757 | parentarray_chunk(i,6,1,2),parentarray_chunk(i,6,2,2),.TRUE.,nb,ndir) |
---|
758 | |
---|
759 | ! call procname(tempP%array6, & |
---|
760 | ! parentarray(1,1,2),parentarray(1,2,2), & |
---|
761 | ! parentarray(2,1,2),parentarray(2,2,2), & |
---|
762 | ! parentarray(3,1,2),parentarray(3,2,2), & |
---|
763 | ! parentarray(4,1,2),parentarray(4,2,2), & |
---|
764 | ! parentarray(5,1,2),parentarray(5,2,2), & |
---|
765 | ! parentarray(6,1,2),parentarray(6,2,2),.TRUE.,nb,ndir) |
---|
766 | end select |
---|
767 | ! |
---|
768 | ! |
---|
769 | endif |
---|
770 | enddo |
---|
771 | call Agrif_ParentGrid_to_ChildGrid() |
---|
772 | nullify(Agrif_CurChildgrid) |
---|
773 | |
---|
774 | #if defined AGRIF_MPI |
---|
775 | if (.not.find_list_interp) then |
---|
776 | ! |
---|
777 | tab3(:,1) = indminglob2(:) |
---|
778 | tab3(:,2) = indmaxglob2(:) |
---|
779 | tab3(:,3) = indmin(:) |
---|
780 | tab3(:,4) = indmax(:) |
---|
781 | tab5(:,1) = indminglob3(:) |
---|
782 | tab5(:,2) = indmaxglob3(:) |
---|
783 | if (agrif_debug_interp) then |
---|
784 | print *,'********************' |
---|
785 | print *,'MPI VARIABLES' |
---|
786 | print *,'INDMINGLOB2' |
---|
787 | do i=1,nbdim |
---|
788 | print *,'Direction ',i,indminglob2(i),indmaxglob2(i) |
---|
789 | enddo |
---|
790 | print *,'INDMIN' |
---|
791 | do i=1,nbdim |
---|
792 | print *,'Direction ',i,indmin(i),indmax(i) |
---|
793 | enddo |
---|
794 | print *,'INDMINGLOB3' |
---|
795 | do i=1,nbdim |
---|
796 | print *,'Direction ',i,indminglob3(i),indmaxglob3(i) |
---|
797 | enddo |
---|
798 | endif |
---|
799 | ! |
---|
800 | call MPI_ALLGATHER(tab3,4*nbdim,MPI_INTEGER,tab4,4*nbdim,MPI_INTEGER,Agrif_mpi_comm,code) |
---|
801 | call MPI_ALLGATHER(tab5,2*nbdim,MPI_INTEGER,tab6,2*nbdim,MPI_INTEGER,Agrif_mpi_comm,code) |
---|
802 | if (.not.associated(tempPextend)) allocate(tempPextend) |
---|
803 | |
---|
804 | do k=0,Agrif_Nbprocs-1 |
---|
805 | do j=1,4 |
---|
806 | do i=1,nbdim |
---|
807 | tab4t(i,k,j) = tab4(i,j,k) |
---|
808 | enddo |
---|
809 | enddo |
---|
810 | enddo |
---|
811 | |
---|
812 | do k=0,Agrif_Nbprocs-1 |
---|
813 | do j=1,2 |
---|
814 | do i=1,nbdim |
---|
815 | tab5t(i,k,j) = tab6(i,j,k) |
---|
816 | enddo |
---|
817 | enddo |
---|
818 | enddo |
---|
819 | |
---|
820 | memberin1(1) = memberin |
---|
821 | call MPI_ALLGATHER(memberin1,1,MPI_LOGICAL,memberinall,1,MPI_LOGICAL,Agrif_mpi_comm,code) |
---|
822 | |
---|
823 | call Get_External_Data_first(tab4t(:,:,1),tab4t(:,:,2), & |
---|
824 | tab4t(:,:,3),tab4t(:,:,4), & |
---|
825 | nbdim,memberinall, coords, & |
---|
826 | sendtoproc1,recvfromproc1, & |
---|
827 | tab4t(:,:,5),tab4t(:,:,6), & |
---|
828 | tab4t(:,:,7),tab4t(:,:,8), & |
---|
829 | tab5t(:,:,1),tab5t(:,:,2)) |
---|
830 | endif |
---|
831 | |
---|
832 | call ExchangeSameLevel(sendtoproc1,recvfromproc1,nbdim, & |
---|
833 | tab4t(:,:,3),tab4t(:,:,4),tab4t(:,:,5),tab4t(:,:,6), & |
---|
834 | tab4t(:,:,7),tab4t(:,:,8),memberin,tempP,tempPextend) |
---|
835 | #else |
---|
836 | tempPextend => tempP |
---|
837 | #endif |
---|
838 | |
---|
839 | if (.not.find_list_interp) then |
---|
840 | call Agrif_Addto_list_interp( & |
---|
841 | child%list_interp,pttab,petab, & |
---|
842 | pttab_Child,pttab_Parent,indmin,indmax, & |
---|
843 | indminglob,indmaxglob, & |
---|
844 | pttruetab,cetruetab, & |
---|
845 | memberin,nbdim, & |
---|
846 | parentarray_chunk,parentarray_chunk_decal,decal_chunks,& |
---|
847 | correction_required,member_chuncks,nb_chunks & |
---|
848 | #if defined AGRIF_MPI |
---|
849 | ,indminglob2,indmaxglob2, & |
---|
850 | parentarray, & |
---|
851 | member, & |
---|
852 | tab4t,memberinall,sendtoproc1,recvfromproc1 & |
---|
853 | #endif |
---|
854 | ) |
---|
855 | endif |
---|
856 | ! |
---|
857 | if (memberin) then |
---|
858 | ! |
---|
859 | if (.not.associated(tempC)) allocate(tempC) |
---|
860 | ! |
---|
861 | call Agrif_array_allocate(tempC,pttruetab,cetruetab,nbdim) |
---|
862 | ! |
---|
863 | ! Special values on the parent grid |
---|
864 | if (Agrif_UseSpecialValue) then |
---|
865 | ! |
---|
866 | noraftab(1:nbdim) = child % root_var % interptab(1:nbdim) == 'N' |
---|
867 | ! |
---|
868 | if (.not.associated(parentvalues)) allocate(parentvalues) |
---|
869 | ! |
---|
870 | call Agrif_array_allocate(parentvalues,indmin,indmax,nbdim) |
---|
871 | call Agrif_var_full_copy_array(parentvalues,tempPextend,nbdim) |
---|
872 | ! |
---|
873 | call Agrif_CheckMasknD(tempPextend,parentvalues, & |
---|
874 | indmin(1:nbdim),indmax(1:nbdim), & |
---|
875 | indmin(1:nbdim),indmax(1:nbdim), & |
---|
876 | noraftab(1:nbdim),nbdim) |
---|
877 | ! |
---|
878 | call Agrif_array_deallocate(parentvalues,nbdim) |
---|
879 | ! |
---|
880 | endif |
---|
881 | ! |
---|
882 | ! Interpolation of the current grid |
---|
883 | ! |
---|
884 | if ( memberin ) then |
---|
885 | select case(nbdim) |
---|
886 | case(1) |
---|
887 | call Agrif_Interp_1D_recursive( type_interp(1), & |
---|
888 | tempPextend%array1, & |
---|
889 | tempC%array1, & |
---|
890 | indmin(1), indmax(1), & |
---|
891 | pttruetab(1), cetruetab(1), & |
---|
892 | s_Child_temp(1), s_Parent_temp(1), & |
---|
893 | ds_Child(1), ds_Parent(1) ) |
---|
894 | case(2) |
---|
895 | call Agrif_Interp_2D_recursive( type_interp(1:2), & |
---|
896 | tempPextend % array2, & |
---|
897 | tempC % array2, & |
---|
898 | indmin(1:2), indmax(1:2), & |
---|
899 | pttruetab(1:2), cetruetab(1:2), & |
---|
900 | s_Child_temp(1:2), s_Parent_temp(1:2), & |
---|
901 | ds_Child(1:2), ds_Parent(1:2) ) |
---|
902 | case(3) |
---|
903 | if (agrif_debug_interp) then |
---|
904 | print *,'APRES ECHANGE' |
---|
905 | print *,'nombre de chunks ',nb_chunks |
---|
906 | print *,'indmin = ',indmin |
---|
907 | print *,'indmax = ',indmax |
---|
908 | do i=1,nb_chunks |
---|
909 | print *,'CHUNK Number : ',i |
---|
910 | print *,'MEMBER = ',member_chuncks(i) |
---|
911 | print *,'Correction = ',correction_required(i) |
---|
912 | enddo |
---|
913 | endif |
---|
914 | |
---|
915 | if (agrif_debug_interp) then |
---|
916 | ! if ((nb==1).AND.(ndir==1)) then |
---|
917 | print *,'valeur parent = ' |
---|
918 | do j=indmin(2),indmax(2) |
---|
919 | do i=indmin(1),indmax(1) |
---|
920 | print *,'par = ',i,j,tempPextend%array3(i,j,1) |
---|
921 | enddo |
---|
922 | enddo |
---|
923 | |
---|
924 | ! endif |
---|
925 | endif |
---|
926 | call Agrif_Interp_3D_recursive( type_interp(1:3), & |
---|
927 | tempPextend % array3, & |
---|
928 | tempC % array3, & |
---|
929 | indmin(1:3), indmax(1:3), & |
---|
930 | pttruetab(1:3), cetruetab(1:3), & |
---|
931 | s_Child_temp(1:3), s_Parent_temp(1:3), & |
---|
932 | ds_Child(1:3), ds_Parent(1:3) ) |
---|
933 | if (agrif_debug_interp) then |
---|
934 | ! if ((nb==1).AND.(ndir==1)) then |
---|
935 | print *,'valeur enfnat = ' |
---|
936 | do j=pttruetab(2),cetruetab(2) |
---|
937 | do i=pttruetab(1),cetruetab(1) |
---|
938 | print *,'par = ',i,j,tempC%array3(i,j,1) |
---|
939 | enddo |
---|
940 | enddo |
---|
941 | |
---|
942 | ! endif |
---|
943 | endif |
---|
944 | case(4) |
---|
945 | call Agrif_Interp_4D_recursive( type_interp(1:4), & |
---|
946 | tempPextend % array4, & |
---|
947 | tempC % array4, & |
---|
948 | indmin(1:4), indmax(1:4), & |
---|
949 | pttruetab(1:4), cetruetab(1:4), & |
---|
950 | s_Child_temp(1:4), s_Parent_temp(1:4), & |
---|
951 | ds_Child(1:4), ds_Parent(1:4) ) |
---|
952 | case(5) |
---|
953 | call Agrif_Interp_5D_recursive( type_interp(1:5), & |
---|
954 | tempPextend % array5, & |
---|
955 | tempC % array5, & |
---|
956 | indmin(1:5), indmax(1:5), & |
---|
957 | pttruetab(1:5), cetruetab(1:5), & |
---|
958 | s_Child_temp(1:5), s_Parent_temp(1:5), & |
---|
959 | ds_Child(1:5), ds_Parent(1:5) ) |
---|
960 | case(6) |
---|
961 | call Agrif_Interp_6D_recursive( type_interp(1:6), & |
---|
962 | tempPextend % array6, & |
---|
963 | tempC % array6, & |
---|
964 | indmin(1:6), indmax(1:6), & |
---|
965 | pttruetab(1:6), cetruetab(1:6), & |
---|
966 | s_Child_temp(1:6), s_Parent_temp(1:6), & |
---|
967 | ds_Child(1:6), ds_Parent(1:6) ) |
---|
968 | end select |
---|
969 | ! |
---|
970 | call Agrif_get_var_bounds_array(child,lowerbound,upperbound,nbdim,parent) |
---|
971 | |
---|
972 | #if defined AGRIF_MPI |
---|
973 | call Agrif_GlobalToLocalBounds(childarray, lowerbound, upperbound, & |
---|
974 | pttruetab, cetruetab, coords, & |
---|
975 | nbdim, Agrif_Procrank, memberout) |
---|
976 | #else |
---|
977 | childarray(:,1,1) = pttruetab |
---|
978 | childarray(:,2,1) = cetruetab |
---|
979 | childarray(:,1,2) = pttruetab |
---|
980 | childarray(:,2,2) = cetruetab |
---|
981 | !cccccccccccccc memberout = .TRUE. |
---|
982 | #endif |
---|
983 | ! |
---|
984 | ! Special values on the child grid |
---|
985 | if (Agrif_UseSpecialValueFineGrid) then |
---|
986 | call GiveAgrif_SpecialValueToTab_mpi( child, tempC, childarray, Agrif_SpecialValueFineGrid,nbdim ) |
---|
987 | endif |
---|
988 | ! |
---|
989 | endif ! ( memberin ) |
---|
990 | ! |
---|
991 | if (torestore) then |
---|
992 | ! |
---|
993 | #if defined AGRIF_MPI |
---|
994 | ! |
---|
995 | SELECT CASE (nbdim) |
---|
996 | CASE (1) |
---|
997 | do i = pttruetab(1),cetruetab(1) |
---|
998 | !hildarrayAModifier if (restore%restore1D(i) == 0) & |
---|
999 | !hildarrayAModifier child%array1(childarray(i,1,2)) = tempC%array1(i) |
---|
1000 | enddo |
---|
1001 | CASE (2) |
---|
1002 | do i = pttruetab(1),cetruetab(1) |
---|
1003 | do j = pttruetab(2),cetruetab(2) |
---|
1004 | !hildarrayAModifier if (restore%restore2D(i,j) == 0) & |
---|
1005 | !hildarrayAModifier child%array2(childarray(i,1,2), & |
---|
1006 | !hildarrayAModifier childarray(j,2,2)) = tempC%array2(i,j) |
---|
1007 | enddo |
---|
1008 | enddo |
---|
1009 | CASE (3) |
---|
1010 | do i = pttruetab(1),cetruetab(1) |
---|
1011 | do j = pttruetab(2),cetruetab(2) |
---|
1012 | do k = pttruetab(3),cetruetab(3) |
---|
1013 | !hildarrayAModifier if (restore%restore3D(i,j,k) == 0) & |
---|
1014 | !hildarrayAModifier child%array3(childarray(i,1,2), & |
---|
1015 | !hildarrayAModifier childarray(j,2,2), & |
---|
1016 | !hildarrayAModifier childarray(k,3,2)) = tempC%array3(i,j,k) |
---|
1017 | enddo |
---|
1018 | enddo |
---|
1019 | enddo |
---|
1020 | CASE (4) |
---|
1021 | do i = pttruetab(1),cetruetab(1) |
---|
1022 | do j = pttruetab(2),cetruetab(2) |
---|
1023 | do k = pttruetab(3),cetruetab(3) |
---|
1024 | do l = pttruetab(4),cetruetab(4) |
---|
1025 | !hildarrayAModifier if (restore%restore4D(i,j,k,l) == 0) & |
---|
1026 | !hildarrayAModifier child%array4(childarray(i,1,2), & |
---|
1027 | !hildarrayAModifier childarray(j,2,2), & |
---|
1028 | !hildarrayAModifier childarray(k,3,2), & |
---|
1029 | !hildarrayAModifier childarray(l,4,2)) = tempC%array4(i,j,k,l) |
---|
1030 | enddo |
---|
1031 | enddo |
---|
1032 | enddo |
---|
1033 | enddo |
---|
1034 | CASE (5) |
---|
1035 | do i = pttruetab(1),cetruetab(1) |
---|
1036 | do j = pttruetab(2),cetruetab(2) |
---|
1037 | do k = pttruetab(3),cetruetab(3) |
---|
1038 | do l = pttruetab(4),cetruetab(4) |
---|
1039 | do m = pttruetab(5),cetruetab(5) |
---|
1040 | !hildarrayAModifier if (restore%restore5D(i,j,k,l,m) == 0) & |
---|
1041 | !hildarrayAModifier child%array5(childarray(i,1,2), & |
---|
1042 | !hildarrayAModifier childarray(j,2,2), & |
---|
1043 | !hildarrayAModifier childarray(k,3,2), & |
---|
1044 | !hildarrayAModifier childarray(l,4,2), & |
---|
1045 | !hildarrayAModifier childarray(m,5,2)) = tempC%array5(i,j,k,l,m) |
---|
1046 | enddo |
---|
1047 | enddo |
---|
1048 | enddo |
---|
1049 | enddo |
---|
1050 | enddo |
---|
1051 | CASE (6) |
---|
1052 | do i = pttruetab(1),cetruetab(1) |
---|
1053 | do j = pttruetab(2),cetruetab(2) |
---|
1054 | do k = pttruetab(3),cetruetab(3) |
---|
1055 | do l = pttruetab(4),cetruetab(4) |
---|
1056 | do m = pttruetab(5),cetruetab(5) |
---|
1057 | do n = pttruetab(6),cetruetab(6) |
---|
1058 | !hildarrayAModifier if (restore%restore6D(i,j,k,l,m,n) == 0) & |
---|
1059 | !hildarrayAModifier child%array6(childarray(i,1,2), & |
---|
1060 | !hildarrayAModifier childarray(j,2,2), & |
---|
1061 | !hildarrayAModifier childarray(k,3,2), & |
---|
1062 | !hildarrayAModifier childarray(l,4,2), & |
---|
1063 | !hildarrayAModifier childarray(m,5,2), & |
---|
1064 | !hildarrayAModifier childarray(n,6,2)) = tempC%array6(i,j,k,l,m,n) |
---|
1065 | enddo |
---|
1066 | enddo |
---|
1067 | enddo |
---|
1068 | enddo |
---|
1069 | enddo |
---|
1070 | enddo |
---|
1071 | END SELECT |
---|
1072 | ! |
---|
1073 | #else |
---|
1074 | select case (nbdim) |
---|
1075 | case (1) |
---|
1076 | do i = pttruetab(1),cetruetab(1) |
---|
1077 | if (restore%restore1D(i) == 0) & |
---|
1078 | parray1(i) = tempC % array1(i) |
---|
1079 | enddo |
---|
1080 | case (2) |
---|
1081 | do j = pttruetab(2),cetruetab(2) |
---|
1082 | do i = pttruetab(1),cetruetab(1) |
---|
1083 | if (restore%restore2D(i,j) == 0) & |
---|
1084 | parray2(i,j) = tempC % array2(i,j) |
---|
1085 | enddo |
---|
1086 | enddo |
---|
1087 | case (3) |
---|
1088 | do k = pttruetab(3),cetruetab(3) |
---|
1089 | do j = pttruetab(2),cetruetab(2) |
---|
1090 | do i = pttruetab(1),cetruetab(1) |
---|
1091 | if (restore%restore3D(i,j,k) == 0) & |
---|
1092 | parray3(i,j,k) = tempC % array3(i,j,k) |
---|
1093 | enddo |
---|
1094 | enddo |
---|
1095 | enddo |
---|
1096 | case (4) |
---|
1097 | do l = pttruetab(4),cetruetab(4) |
---|
1098 | do k = pttruetab(3),cetruetab(3) |
---|
1099 | do j = pttruetab(2),cetruetab(2) |
---|
1100 | do i = pttruetab(1),cetruetab(1) |
---|
1101 | if (restore%restore4D(i,j,k,l) == 0) & |
---|
1102 | parray4(i,j,k,l) = tempC % array4(i,j,k,l) |
---|
1103 | enddo |
---|
1104 | enddo |
---|
1105 | enddo |
---|
1106 | enddo |
---|
1107 | case (5) |
---|
1108 | do m = pttruetab(5),cetruetab(5) |
---|
1109 | do l = pttruetab(4),cetruetab(4) |
---|
1110 | do k = pttruetab(3),cetruetab(3) |
---|
1111 | do j = pttruetab(2),cetruetab(2) |
---|
1112 | do i = pttruetab(1),cetruetab(1) |
---|
1113 | if (restore%restore5D(i,j,k,l,m) == 0) & |
---|
1114 | parray5(i,j,k,l,m) = tempC % array5(i,j,k,l,m) |
---|
1115 | enddo |
---|
1116 | enddo |
---|
1117 | enddo |
---|
1118 | enddo |
---|
1119 | enddo |
---|
1120 | case (6) |
---|
1121 | do n = pttruetab(6),cetruetab(6) |
---|
1122 | do m = pttruetab(5),cetruetab(5) |
---|
1123 | do l = pttruetab(4),cetruetab(4) |
---|
1124 | do k = pttruetab(3),cetruetab(3) |
---|
1125 | do j = pttruetab(2),cetruetab(2) |
---|
1126 | do i = pttruetab(1),cetruetab(1) |
---|
1127 | if (restore%restore6D(i,j,k,l,m,n) == 0) & |
---|
1128 | parray6(i,j,k,l,m,n) = tempC % array6(i,j,k,l,m,n) |
---|
1129 | enddo |
---|
1130 | enddo |
---|
1131 | enddo |
---|
1132 | enddo |
---|
1133 | enddo |
---|
1134 | enddo |
---|
1135 | end select |
---|
1136 | ! |
---|
1137 | #endif |
---|
1138 | ! |
---|
1139 | else ! .not.to_restore |
---|
1140 | ! |
---|
1141 | if (memberin) then |
---|
1142 | ! |
---|
1143 | if ( .not.in_bc ) then |
---|
1144 | select case(nbdim) |
---|
1145 | case(1) |
---|
1146 | call procname(tempC % array1( & |
---|
1147 | childarray(1,1,1):childarray(1,2,1)), & |
---|
1148 | childarray(1,1,2),childarray(1,2,2),.FALSE.,nb,ndir) |
---|
1149 | case(2) |
---|
1150 | call procname( & |
---|
1151 | tempC % array2( & |
---|
1152 | childarray(1,1,1):childarray(1,2,1), & |
---|
1153 | childarray(2,1,1):childarray(2,2,1)), & |
---|
1154 | childarray(1,1,2),childarray(1,2,2), & |
---|
1155 | childarray(2,1,2),childarray(2,2,2),.FALSE.,nb,ndir) |
---|
1156 | case(3) |
---|
1157 | call procname( & |
---|
1158 | tempC % array3( & |
---|
1159 | childarray(1,1,1):childarray(1,2,1), & |
---|
1160 | childarray(2,1,1):childarray(2,2,1), & |
---|
1161 | childarray(3,1,1):childarray(3,2,1)), & |
---|
1162 | childarray(1,1,2),childarray(1,2,2), & |
---|
1163 | childarray(2,1,2),childarray(2,2,2), & |
---|
1164 | childarray(3,1,2),childarray(3,2,2),.FALSE.,nb,ndir) |
---|
1165 | case(4) |
---|
1166 | call procname( & |
---|
1167 | tempC % array4( & |
---|
1168 | childarray(1,1,1):childarray(1,2,1), & |
---|
1169 | childarray(2,1,1):childarray(2,2,1), & |
---|
1170 | childarray(3,1,1):childarray(3,2,1), & |
---|
1171 | childarray(4,1,1):childarray(4,2,1)), & |
---|
1172 | childarray(1,1,2),childarray(1,2,2), & |
---|
1173 | childarray(2,1,2),childarray(2,2,2), & |
---|
1174 | childarray(3,1,2),childarray(3,2,2), & |
---|
1175 | childarray(4,1,2),childarray(4,2,2),.FALSE.,nb,ndir) |
---|
1176 | case(5) |
---|
1177 | call procname( & |
---|
1178 | tempC % array5( & |
---|
1179 | childarray(1,1,1):childarray(1,2,1), & |
---|
1180 | childarray(2,1,1):childarray(2,2,1), & |
---|
1181 | childarray(3,1,1):childarray(3,2,1), & |
---|
1182 | childarray(4,1,1):childarray(4,2,1), & |
---|
1183 | childarray(5,1,1):childarray(5,2,1)), & |
---|
1184 | childarray(1,1,2),childarray(1,2,2), & |
---|
1185 | childarray(2,1,2),childarray(2,2,2), & |
---|
1186 | childarray(3,1,2),childarray(3,2,2), & |
---|
1187 | childarray(4,1,2),childarray(4,2,2), & |
---|
1188 | childarray(5,1,2),childarray(5,2,2),.FALSE.,nb,ndir) |
---|
1189 | case(6) |
---|
1190 | call procname( & |
---|
1191 | tempC % array6( & |
---|
1192 | childarray(1,1,1):childarray(1,2,1), & |
---|
1193 | childarray(2,1,1):childarray(2,2,1), & |
---|
1194 | childarray(3,1,1):childarray(3,2,1), & |
---|
1195 | childarray(4,1,1):childarray(4,2,1), & |
---|
1196 | childarray(5,1,1):childarray(5,2,1), & |
---|
1197 | childarray(6,1,1):childarray(6,2,1)), & |
---|
1198 | childarray(1,1,2),childarray(1,2,2), & |
---|
1199 | childarray(2,1,2),childarray(2,2,2), & |
---|
1200 | childarray(3,1,2),childarray(3,2,2), & |
---|
1201 | childarray(4,1,2),childarray(4,2,2), & |
---|
1202 | childarray(5,1,2),childarray(5,2,2), & |
---|
1203 | childarray(6,1,2),childarray(6,2,2),.FALSE.,nb,ndir) |
---|
1204 | end select |
---|
1205 | else ! we are in_bc |
---|
1206 | select case (nbdim) |
---|
1207 | case (1) |
---|
1208 | parray1(childarray(1,1,2):childarray(1,2,2)) = & |
---|
1209 | tempC%array1(childarray(1,1,1):childarray(1,2,1)) |
---|
1210 | case (2) |
---|
1211 | parray2(childarray(1,1,2):childarray(1,2,2), & |
---|
1212 | childarray(2,1,2):childarray(2,2,2)) = & |
---|
1213 | tempC%array2(childarray(1,1,1):childarray(1,2,1), & |
---|
1214 | childarray(2,1,1):childarray(2,2,1)) |
---|
1215 | case (3) |
---|
1216 | parray3(childarray(1,1,2):childarray(1,2,2), & |
---|
1217 | childarray(2,1,2):childarray(2,2,2), & |
---|
1218 | childarray(3,1,2):childarray(3,2,2)) = & |
---|
1219 | tempC%array3(childarray(1,1,1):childarray(1,2,1), & |
---|
1220 | childarray(2,1,1):childarray(2,2,1), & |
---|
1221 | childarray(3,1,1):childarray(3,2,1)) |
---|
1222 | if (agrif_debug_interp) then |
---|
1223 | if ((nb==1).AND.(ndir==1)) then |
---|
1224 | print *,'valeur enfnat2 = ' |
---|
1225 | do j=childarray(2,1,2),childarray(2,2,2) |
---|
1226 | do i=childarray(1,1,2),childarray(1,2,2) |
---|
1227 | print *,'par = ',i,j,parray3(i,j,1) |
---|
1228 | enddo |
---|
1229 | enddo |
---|
1230 | |
---|
1231 | endif |
---|
1232 | endif |
---|
1233 | case (4) |
---|
1234 | parray4(childarray(1,1,2):childarray(1,2,2), & |
---|
1235 | childarray(2,1,2):childarray(2,2,2), & |
---|
1236 | childarray(3,1,2):childarray(3,2,2), & |
---|
1237 | childarray(4,1,2):childarray(4,2,2)) = & |
---|
1238 | tempC%array4(childarray(1,1,1):childarray(1,2,1), & |
---|
1239 | childarray(2,1,1):childarray(2,2,1), & |
---|
1240 | childarray(3,1,1):childarray(3,2,1), & |
---|
1241 | childarray(4,1,1):childarray(4,2,1)) |
---|
1242 | case (5) |
---|
1243 | parray5(childarray(1,1,2):childarray(1,2,2), & |
---|
1244 | childarray(2,1,2):childarray(2,2,2), & |
---|
1245 | childarray(3,1,2):childarray(3,2,2), & |
---|
1246 | childarray(4,1,2):childarray(4,2,2), & |
---|
1247 | childarray(5,1,2):childarray(5,2,2)) = & |
---|
1248 | tempC%array5(childarray(1,1,1):childarray(1,2,1), & |
---|
1249 | childarray(2,1,1):childarray(2,2,1), & |
---|
1250 | childarray(3,1,1):childarray(3,2,1), & |
---|
1251 | childarray(4,1,1):childarray(4,2,1), & |
---|
1252 | childarray(5,1,1):childarray(5,2,1)) |
---|
1253 | case (6) |
---|
1254 | parray6(childarray(1,1,2):childarray(1,2,2), & |
---|
1255 | childarray(2,1,2):childarray(2,2,2), & |
---|
1256 | childarray(3,1,2):childarray(3,2,2), & |
---|
1257 | childarray(4,1,2):childarray(4,2,2), & |
---|
1258 | childarray(5,1,2):childarray(5,2,2), & |
---|
1259 | childarray(6,1,2):childarray(6,2,2)) = & |
---|
1260 | tempC%array6(childarray(1,1,1):childarray(1,2,1), & |
---|
1261 | childarray(2,1,1):childarray(2,2,1), & |
---|
1262 | childarray(3,1,1):childarray(3,2,1), & |
---|
1263 | childarray(4,1,1):childarray(4,2,1), & |
---|
1264 | childarray(5,1,1):childarray(5,2,1), & |
---|
1265 | childarray(6,1,1):childarray(6,2,1)) |
---|
1266 | end select |
---|
1267 | endif ! < (.not.in_bc) |
---|
1268 | endif ! < memberin |
---|
1269 | ! |
---|
1270 | endif |
---|
1271 | |
---|
1272 | call Agrif_array_deallocate(tempPextend,nbdim) |
---|
1273 | call Agrif_array_deallocate(tempC,nbdim) |
---|
1274 | |
---|
1275 | endif |
---|
1276 | ! |
---|
1277 | ! Deallocations |
---|
1278 | #if defined AGRIF_MPI |
---|
1279 | if (member) then |
---|
1280 | if (agrif_debug_interp) then |
---|
1281 | print *,'ALLCOATED 0 = ',allocated(tempP%array3),size(tempP%array3) |
---|
1282 | endif |
---|
1283 | call Agrif_array_deallocate(tempP,nbdim) |
---|
1284 | endif |
---|
1285 | #endif |
---|
1286 | !--------------------------------------------------------------------------------------------------- |
---|
1287 | end subroutine Agrif_InterpnD |
---|
1288 | !=================================================================================================== |
---|
1289 | ! |
---|
1290 | !=================================================================================================== |
---|
1291 | ! subroutine Agrif_Parentbounds |
---|
1292 | ! |
---|
1293 | !> Calculates the bounds of the parent grid for the interpolation of the child grid |
---|
1294 | !--------------------------------------------------------------------------------------------------- |
---|
1295 | subroutine Agrif_Parentbounds ( type_interp, nbdim, indmin, indmax, & |
---|
1296 | s_Parent_temp, s_Child_temp, & |
---|
1297 | s_Child, ds_Child, & |
---|
1298 | s_Parent,ds_Parent, & |
---|
1299 | pttruetab, cetruetab, & |
---|
1300 | pttab_Child, pttab_Parent, posvar, coords ) |
---|
1301 | !--------------------------------------------------------------------------------------------------- |
---|
1302 | INTEGER, DIMENSION(6), intent(in) :: type_interp |
---|
1303 | INTEGER, intent(in) :: nbdim |
---|
1304 | INTEGER, DIMENSION(nbdim), intent(out) :: indmin, indmax |
---|
1305 | REAL, DIMENSION(nbdim), intent(out) :: s_Parent_temp, s_child_temp |
---|
1306 | REAL, DIMENSION(nbdim), intent(in) :: s_Child, ds_child |
---|
1307 | REAL, DIMENSION(nbdim), intent(in) :: s_Parent,ds_Parent |
---|
1308 | INTEGER, DIMENSION(nbdim), intent(in) :: pttruetab, cetruetab |
---|
1309 | INTEGER, DIMENSION(nbdim), intent(in) :: pttab_Child, pttab_Parent |
---|
1310 | INTEGER, DIMENSION(nbdim), intent(in) :: posvar |
---|
1311 | INTEGER, DIMENSION(nbdim), intent(in) :: coords |
---|
1312 | ! |
---|
1313 | INTEGER :: i |
---|
1314 | REAL,DIMENSION(nbdim) :: dim_newmin, dim_newmax |
---|
1315 | ! |
---|
1316 | dim_newmin = s_Child + (pttruetab - pttab_Child) * ds_Child |
---|
1317 | dim_newmax = s_Child + (cetruetab - pttab_Child) * ds_Child |
---|
1318 | ! |
---|
1319 | do i = 1,nbdim |
---|
1320 | ! |
---|
1321 | indmin(i) = pttab_Parent(i) + agrif_int((dim_newmin(i)-s_Parent(i))/ds_Parent(i)) |
---|
1322 | indmax(i) = pttab_Parent(i) + agrif_ceiling((dim_newmax(i)-s_Parent(i))/ds_Parent(i)) |
---|
1323 | ! |
---|
1324 | ! Necessary for the Quadratic interpolation |
---|
1325 | ! |
---|
1326 | if ( (pttruetab(i) == cetruetab(i)) .and. (posvar(i) == 1) ) then |
---|
1327 | elseif ( coords(i) == 0 ) then ! (interptab == 'N') |
---|
1328 | elseif ( (type_interp(i) == Agrif_ppm) .or. & |
---|
1329 | (type_interp(i) == Agrif_eno) .or. & |
---|
1330 | (type_interp(i) == Agrif_ppm_lim) .or. & |
---|
1331 | (type_interp(i) == Agrif_weno) ) then |
---|
1332 | indmin(i) = indmin(i) - 2 |
---|
1333 | indmax(i) = indmax(i) + 2 |
---|
1334 | |
---|
1335 | if (Agrif_UseSpecialValue) then |
---|
1336 | indmin(i) = indmin(i)-MaxSearch |
---|
1337 | indmax(i) = indmax(i)+MaxSearch |
---|
1338 | endif |
---|
1339 | |
---|
1340 | elseif ( (type_interp(i) /= Agrif_constant) .and. & |
---|
1341 | (type_interp(i) /= Agrif_linear) ) then |
---|
1342 | indmin(i) = indmin(i) - 1 |
---|
1343 | indmax(i) = indmax(i) + 1 |
---|
1344 | |
---|
1345 | if (Agrif_UseSpecialValue) then |
---|
1346 | indmin(i) = indmin(i)-MaxSearch |
---|
1347 | indmax(i) = indmax(i)+MaxSearch |
---|
1348 | endif |
---|
1349 | |
---|
1350 | elseif ( (type_interp(i) == Agrif_constant) .or. & |
---|
1351 | (type_interp(i) == Agrif_linear) ) then |
---|
1352 | if (Agrif_UseSpecialValue) then |
---|
1353 | indmin(i) = indmin(i)-MaxSearch |
---|
1354 | indmax(i) = indmax(i)+MaxSearch |
---|
1355 | endif |
---|
1356 | |
---|
1357 | endif |
---|
1358 | ! |
---|
1359 | enddo |
---|
1360 | ! |
---|
1361 | s_Parent_temp = s_Parent + (indmin - pttab_Parent) * ds_Parent |
---|
1362 | s_Child_temp = s_Child + (pttruetab - pttab_Child) * ds_Child |
---|
1363 | !--------------------------------------------------------------------------------------------------- |
---|
1364 | end subroutine Agrif_Parentbounds |
---|
1365 | !=================================================================================================== |
---|
1366 | ! |
---|
1367 | !=================================================================================================== |
---|
1368 | ! subroutine Agrif_Interp_1D_Recursive |
---|
1369 | ! |
---|
1370 | !> Subroutine for the interpolation of a 1D grid variable. |
---|
1371 | !> It calls Agrif_InterpBase. |
---|
1372 | !--------------------------------------------------------------------------------------------------- |
---|
1373 | subroutine Agrif_Interp_1D_recursive ( type_interp, tabin, tabout, & |
---|
1374 | indmin, indmax, & |
---|
1375 | pttab_child, petab_child, & |
---|
1376 | s_child, s_parent, & |
---|
1377 | ds_child, ds_parent ) |
---|
1378 | !--------------------------------------------------------------------------------------------------- |
---|
1379 | integer, intent(in) :: type_interp |
---|
1380 | integer, intent(in) :: indmin, indmax |
---|
1381 | integer, intent(in) :: pttab_child, petab_child |
---|
1382 | real, intent(in) :: s_child, s_parent |
---|
1383 | real, intent(in) :: ds_child, ds_parent |
---|
1384 | real, dimension( & |
---|
1385 | indmin:indmax & |
---|
1386 | ), intent(in) :: tabin |
---|
1387 | real, dimension( & |
---|
1388 | pttab_child:petab_child & |
---|
1389 | ), intent(out) :: tabout |
---|
1390 | !--------------------------------------------------------------------------------------------------- |
---|
1391 | call Agrif_InterpBase(type_interp, & |
---|
1392 | tabin(indmin:indmax), & |
---|
1393 | tabout(pttab_child:petab_child), & |
---|
1394 | indmin, indmax, & |
---|
1395 | pttab_child, petab_child, & |
---|
1396 | s_parent, s_child, & |
---|
1397 | ds_parent, ds_child) |
---|
1398 | !--------------------------------------------------------------------------------------------------- |
---|
1399 | end subroutine Agrif_Interp_1D_recursive |
---|
1400 | !=================================================================================================== |
---|
1401 | ! |
---|
1402 | !=================================================================================================== |
---|
1403 | ! subroutine Agrif_Interp_2D_Recursive |
---|
1404 | ! |
---|
1405 | !> Subroutine for the interpolation of a 2D grid variable. |
---|
1406 | !> It calls Agrif_Interp_1D_recursive and Agrif_InterpBase. |
---|
1407 | !--------------------------------------------------------------------------------------------------- |
---|
1408 | subroutine Agrif_Interp_2D_recursive ( type_interp, tabin, tabout, & |
---|
1409 | indmin, indmax, & |
---|
1410 | pttab_child, petab_child, & |
---|
1411 | s_child, s_parent, & |
---|
1412 | ds_child, ds_parent ) |
---|
1413 | !--------------------------------------------------------------------------------------------------- |
---|
1414 | integer, dimension(2), intent(in) :: type_interp |
---|
1415 | integer, dimension(2), intent(in) :: indmin, indmax |
---|
1416 | integer, dimension(2), intent(in) :: pttab_child, petab_child |
---|
1417 | real, dimension(2), intent(in) :: s_child, s_parent |
---|
1418 | real, dimension(2), intent(in) :: ds_child, ds_parent |
---|
1419 | real, dimension( & |
---|
1420 | indmin(1):indmax(1), & |
---|
1421 | indmin(2):indmax(2)), intent(in) :: tabin |
---|
1422 | real, dimension( & |
---|
1423 | pttab_child(1):petab_child(1), & |
---|
1424 | pttab_child(2):petab_child(2)), intent(out) :: tabout |
---|
1425 | !--------------------------------------------------------------------------------------------------- |
---|
1426 | real, dimension( & |
---|
1427 | pttab_child(1):petab_child(1), & |
---|
1428 | indmin(2):indmax(2)) :: tabtemp |
---|
1429 | real, dimension( & |
---|
1430 | pttab_child(2):petab_child(2), & |
---|
1431 | pttab_child(1):petab_child(1)) :: tabout_trsp |
---|
1432 | real, dimension( & |
---|
1433 | indmin(2):indmax(2), & |
---|
1434 | pttab_child(1):petab_child(1)) :: tabtemp_trsp |
---|
1435 | integer :: i, j, coeffraf |
---|
1436 | !--------------------------------------------------------------------------------------------------- |
---|
1437 | ! |
---|
1438 | coeffraf = nint ( ds_parent(1) / ds_child(1) ) |
---|
1439 | ! |
---|
1440 | if ((type_interp(1) == Agrif_Linear) .and. (coeffraf /= 1)) then |
---|
1441 | !---CDIR NEXPAND |
---|
1442 | if(.NOT. precomputedone(1)) & |
---|
1443 | call Linear1dPrecompute2d( & |
---|
1444 | indmax(2)-indmin(2)+1, & |
---|
1445 | indmax(1)-indmin(1)+1, & |
---|
1446 | petab_child(1)-pttab_child(1)+1, & |
---|
1447 | s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
1448 | !---CDIR NEXPAND |
---|
1449 | call Linear1dAfterCompute(tabin,tabtemp,size(tabin),size(tabtemp),1) |
---|
1450 | ! |
---|
1451 | elseif ((type_interp(1) == Agrif_PPM) .and. (coeffraf /= 1)) then |
---|
1452 | !---CDIR NEXPAND |
---|
1453 | if(.NOT. precomputedone(1)) & |
---|
1454 | call PPM1dPrecompute2d( & |
---|
1455 | indmax(2)-indmin(2)+1, & |
---|
1456 | indmax(1)-indmin(1)+1, & |
---|
1457 | petab_child(1)-pttab_child(1)+1, & |
---|
1458 | s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
1459 | !---CDIR NEXPAND |
---|
1460 | call PPM1dAfterCompute(tabin,tabtemp,size(tabin),size(tabtemp),1) |
---|
1461 | else |
---|
1462 | do j = indmin(2),indmax(2) |
---|
1463 | ! |
---|
1464 | !---CDIR NEXPAND |
---|
1465 | call Agrif_Interp_1D_recursive(type_interp(1), & |
---|
1466 | tabin(indmin(1):indmax(1),j), & |
---|
1467 | tabtemp(pttab_child(1):petab_child(1),j), & |
---|
1468 | indmin(1),indmax(1), & |
---|
1469 | pttab_child(1),petab_child(1), & |
---|
1470 | s_child(1), s_parent(1), & |
---|
1471 | ds_child(1),ds_parent(1)) |
---|
1472 | ! |
---|
1473 | enddo |
---|
1474 | endif |
---|
1475 | |
---|
1476 | coeffraf = nint(ds_parent(2)/ds_child(2)) |
---|
1477 | tabtemp_trsp = TRANSPOSE(tabtemp) |
---|
1478 | |
---|
1479 | if ((type_interp(2) == Agrif_Linear) .and. (coeffraf /= 1)) then |
---|
1480 | !---CDIR NEXPAND |
---|
1481 | if(.NOT. precomputedone(2)) & |
---|
1482 | call Linear1dPrecompute2d( & |
---|
1483 | petab_child(1)-pttab_child(1)+1, & |
---|
1484 | indmax(2)-indmin(2)+1, & |
---|
1485 | petab_child(2)-pttab_child(2)+1, & |
---|
1486 | s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
1487 | !---CDIR NEXPAND |
---|
1488 | call Linear1dAfterCompute(tabtemp_trsp,tabout_trsp, & |
---|
1489 | size(tabtemp_trsp),size(tabout_trsp),2) |
---|
1490 | |
---|
1491 | elseif ((type_interp(2) == Agrif_PPM) .and. (coeffraf /= 1)) then |
---|
1492 | !---CDIR NEXPAND |
---|
1493 | if(.NOT. precomputedone(2)) & |
---|
1494 | call PPM1dPrecompute2d( & |
---|
1495 | petab_child(1)-pttab_child(1)+1, & |
---|
1496 | indmax(2)-indmin(2)+1, & |
---|
1497 | petab_child(2)-pttab_child(2)+1, & |
---|
1498 | s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
1499 | !---CDIR NEXPAND |
---|
1500 | call PPM1dAfterCompute(tabtemp_trsp, tabout_trsp, & |
---|
1501 | size(tabtemp_trsp), size(tabout_trsp), 2) |
---|
1502 | else |
---|
1503 | do i = pttab_child(1), petab_child(1) |
---|
1504 | ! |
---|
1505 | !---CDIR NEXPAND |
---|
1506 | call Agrif_InterpBase(type_interp(2), & |
---|
1507 | tabtemp_trsp(indmin(2):indmax(2), i), & |
---|
1508 | tabout_trsp(pttab_child(2):petab_child(2), i), & |
---|
1509 | indmin(2), indmax(2), & |
---|
1510 | pttab_child(2), petab_child(2), & |
---|
1511 | s_parent(2), s_child(2), & |
---|
1512 | ds_parent(2), ds_child(2) ) |
---|
1513 | enddo |
---|
1514 | endif |
---|
1515 | ! |
---|
1516 | tabout = TRANSPOSE(tabout_trsp) |
---|
1517 | !--------------------------------------------------------------------------------------------------- |
---|
1518 | end subroutine Agrif_Interp_2D_recursive |
---|
1519 | !=================================================================================================== |
---|
1520 | ! |
---|
1521 | !=================================================================================================== |
---|
1522 | ! subroutine Agrif_Interp_3D_Recursive |
---|
1523 | ! |
---|
1524 | !> Subroutine for the interpolation of a 3D grid variable. |
---|
1525 | !> It calls #Agrif_Interp_2D_recursive and #Agrif_InterpBase. |
---|
1526 | !--------------------------------------------------------------------------------------------------- |
---|
1527 | subroutine Agrif_Interp_3D_recursive ( type_interp, tabin, tabout, & |
---|
1528 | indmin, indmax, & |
---|
1529 | pttab_child, petab_child, & |
---|
1530 | s_child, s_parent, & |
---|
1531 | ds_child, ds_parent ) |
---|
1532 | !--------------------------------------------------------------------------------------------------- |
---|
1533 | integer, dimension(3), intent(in) :: type_interp |
---|
1534 | integer, dimension(3), intent(in) :: indmin, indmax |
---|
1535 | integer, dimension(3), intent(in) :: pttab_child, petab_child |
---|
1536 | real, dimension(3), intent(in) :: s_child, s_parent |
---|
1537 | real, dimension(3), intent(in) :: ds_child, ds_parent |
---|
1538 | real, dimension( & |
---|
1539 | indmin(1):indmax(1), & |
---|
1540 | indmin(2):indmax(2), & |
---|
1541 | indmin(3):indmax(3)), intent(in) :: tabin |
---|
1542 | real, dimension( & |
---|
1543 | pttab_child(1):petab_child(1), & |
---|
1544 | pttab_child(2):petab_child(2), & |
---|
1545 | pttab_child(3):petab_child(3)), intent(out) :: tabout |
---|
1546 | !--------------------------------------------------------------------------------------------------- |
---|
1547 | real, dimension( & |
---|
1548 | pttab_child(1):petab_child(1), & |
---|
1549 | pttab_child(2):petab_child(2), & |
---|
1550 | indmin(3):indmax(3)) :: tabtemp |
---|
1551 | integer :: i, j, k, coeffraf |
---|
1552 | integer :: locind_child_left, kdeb |
---|
1553 | ! |
---|
1554 | coeffraf = nint ( ds_parent(1) / ds_child(1) ) |
---|
1555 | if ( (type_interp(1) == Agrif_Linear) .and. (coeffraf/=1) ) then |
---|
1556 | call Linear1dPrecompute2d(indmax(2)-indmin(2)+1, & |
---|
1557 | indmax(1)-indmin(1)+1, & |
---|
1558 | petab_child(1)-pttab_child(1)+1, & |
---|
1559 | s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
1560 | precomputedone(1) = .TRUE. |
---|
1561 | elseif ( (type_interp(1) == Agrif_PPM) .and. (coeffraf/=1) ) then |
---|
1562 | call PPM1dPrecompute2d(indmax(2)-indmin(2)+1, & |
---|
1563 | indmax(1)-indmin(1)+1, & |
---|
1564 | petab_child(1)-pttab_child(1)+1, & |
---|
1565 | s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
1566 | precomputedone(1) = .TRUE. |
---|
1567 | endif |
---|
1568 | |
---|
1569 | coeffraf = nint ( ds_parent(2) / ds_child(2) ) |
---|
1570 | if ( (type_interp(2) == Agrif_Linear) .and. (coeffraf/=1) ) then |
---|
1571 | call Linear1dPrecompute2d(petab_child(1)-pttab_child(1)+1, & |
---|
1572 | indmax(2)-indmin(2)+1, & |
---|
1573 | petab_child(2)-pttab_child(2)+1, & |
---|
1574 | s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
1575 | precomputedone(2) = .TRUE. |
---|
1576 | elseif ( (type_interp(2) == Agrif_PPM) .and. (coeffraf/=1) ) then |
---|
1577 | call PPM1dPrecompute2d(petab_child(1)-pttab_child(1)+1, & |
---|
1578 | indmax(2)-indmin(2)+1, & |
---|
1579 | petab_child(2)-pttab_child(2)+1, & |
---|
1580 | s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
1581 | precomputedone(2) = .TRUE. |
---|
1582 | endif |
---|
1583 | ! |
---|
1584 | do k = indmin(3), indmax(3) |
---|
1585 | call Agrif_Interp_2D_recursive(type_interp(1:2), & |
---|
1586 | tabin(indmin(1):indmax(1), & |
---|
1587 | indmin(2):indmax(2), k), & |
---|
1588 | tabtemp(pttab_child(1):petab_child(1), & |
---|
1589 | pttab_child(2):petab_child(2), k), & |
---|
1590 | indmin(1:2), indmax(1:2), & |
---|
1591 | pttab_child(1:2), petab_child(1:2), & |
---|
1592 | s_child(1:2), s_parent(1:2), & |
---|
1593 | ds_child(1:2), ds_parent(1:2) ) |
---|
1594 | enddo |
---|
1595 | ! |
---|
1596 | precomputedone(1) = .FALSE. |
---|
1597 | precomputedone(2) = .FALSE. |
---|
1598 | coeffraf = nint(ds_parent(3)/ds_child(3)) |
---|
1599 | ! |
---|
1600 | if ( coeffraf == 1 ) then |
---|
1601 | locind_child_left = 1 + agrif_int((s_child(3)-s_parent(3))/ds_parent(3)) |
---|
1602 | kdeb = indmin(3)+locind_child_left-2 |
---|
1603 | do k = pttab_child(3),petab_child(3) |
---|
1604 | kdeb = kdeb + 1 |
---|
1605 | do j = pttab_child(2), petab_child(2) |
---|
1606 | do i = pttab_child(1), petab_child(1) |
---|
1607 | tabout(i,j,k) = tabtemp(i,j,kdeb) |
---|
1608 | enddo |
---|
1609 | enddo |
---|
1610 | enddo |
---|
1611 | else |
---|
1612 | do j = pttab_child(2), petab_child(2) |
---|
1613 | do i = pttab_child(1), petab_child(1) |
---|
1614 | call Agrif_InterpBase(type_interp(3), & |
---|
1615 | tabtemp(i,j,indmin(3):indmax(3)), & |
---|
1616 | tabout(i,j,pttab_child(3):petab_child(3)), & |
---|
1617 | indmin(3), indmax(3), & |
---|
1618 | pttab_child(3), petab_child(3), & |
---|
1619 | s_parent(3), s_child(3), & |
---|
1620 | ds_parent(3), ds_child(3) ) |
---|
1621 | enddo |
---|
1622 | enddo |
---|
1623 | endif |
---|
1624 | !--------------------------------------------------------------------------------------------------- |
---|
1625 | end subroutine Agrif_Interp_3D_recursive |
---|
1626 | !=================================================================================================== |
---|
1627 | ! |
---|
1628 | !=================================================================================================== |
---|
1629 | ! subroutine Agrif_Interp_4D_Recursive |
---|
1630 | ! |
---|
1631 | !> Subroutine for the interpolation of a 4D grid variable. |
---|
1632 | !> It calls #Agrif_Interp_3D_recursive and #Agrif_InterpBase. |
---|
1633 | !--------------------------------------------------------------------------------------------------- |
---|
1634 | subroutine Agrif_Interp_4D_recursive ( type_interp, tabin, tabout, & |
---|
1635 | indmin, indmax, & |
---|
1636 | pttab_child, petab_child, & |
---|
1637 | s_child, s_parent, & |
---|
1638 | ds_child, ds_parent ) |
---|
1639 | !--------------------------------------------------------------------------------------------------- |
---|
1640 | integer, dimension(4), intent(in) :: type_interp |
---|
1641 | integer, dimension(4), intent(in) :: indmin, indmax |
---|
1642 | integer, dimension(4), intent(in) :: pttab_child, petab_child |
---|
1643 | real, dimension(4), intent(in) :: s_child, s_parent |
---|
1644 | real, dimension(4), intent(in) :: ds_child, ds_parent |
---|
1645 | real, dimension( & |
---|
1646 | indmin(1):indmax(1), & |
---|
1647 | indmin(2):indmax(2), & |
---|
1648 | indmin(3):indmax(3), & |
---|
1649 | indmin(4):indmax(4)), intent(in) :: tabin |
---|
1650 | real, dimension( & |
---|
1651 | pttab_child(1):petab_child(1), & |
---|
1652 | pttab_child(2):petab_child(2), & |
---|
1653 | pttab_child(3):petab_child(3), & |
---|
1654 | pttab_child(4):petab_child(4)), intent(out) :: tabout |
---|
1655 | !--------------------------------------------------------------------------------------------------- |
---|
1656 | real, dimension( & |
---|
1657 | pttab_child(1):petab_child(1), & |
---|
1658 | pttab_child(2):petab_child(2), & |
---|
1659 | pttab_child(3):petab_child(3), & |
---|
1660 | indmin(4):indmax(4)) :: tabtemp |
---|
1661 | integer :: i, j, k, l |
---|
1662 | ! |
---|
1663 | do l = indmin(4), indmax(4) |
---|
1664 | call Agrif_Interp_3D_recursive(type_interp(1:3), & |
---|
1665 | tabin(indmin(1):indmax(1), & |
---|
1666 | indmin(2):indmax(2), & |
---|
1667 | indmin(3):indmax(3), l), & |
---|
1668 | tabtemp(pttab_child(1):petab_child(1), & |
---|
1669 | pttab_child(2):petab_child(2), & |
---|
1670 | pttab_child(3):petab_child(3), l), & |
---|
1671 | indmin(1:3), indmax(1:3), & |
---|
1672 | pttab_child(1:3), petab_child(1:3), & |
---|
1673 | s_child(1:3), s_parent(1:3), & |
---|
1674 | ds_child(1:3), ds_parent(1:3) ) |
---|
1675 | enddo |
---|
1676 | ! |
---|
1677 | do k = pttab_child(3), petab_child(3) |
---|
1678 | do j = pttab_child(2), petab_child(2) |
---|
1679 | do i = pttab_child(1), petab_child(1) |
---|
1680 | call Agrif_InterpBase(type_interp(4), & |
---|
1681 | tabtemp(i,j,k,indmin(4):indmax(4)), & |
---|
1682 | tabout(i,j,k,pttab_child(4):petab_child(4)), & |
---|
1683 | indmin(4), indmax(4), & |
---|
1684 | pttab_child(4), petab_child(4), & |
---|
1685 | s_parent(4), s_child(4), & |
---|
1686 | ds_parent(4), ds_child(4) ) |
---|
1687 | enddo |
---|
1688 | enddo |
---|
1689 | enddo |
---|
1690 | !--------------------------------------------------------------------------------------------------- |
---|
1691 | end subroutine Agrif_Interp_4D_recursive |
---|
1692 | !=================================================================================================== |
---|
1693 | ! |
---|
1694 | !=================================================================================================== |
---|
1695 | ! subroutine Agrif_Interp_5D_Recursive |
---|
1696 | ! |
---|
1697 | !> Subroutine for the interpolation of a 5D grid variable. |
---|
1698 | !> It calls #Agrif_Interp_4D_recursive and #Agrif_InterpBase. |
---|
1699 | !--------------------------------------------------------------------------------------------------- |
---|
1700 | subroutine Agrif_Interp_5D_recursive ( type_interp, tabin, tabout, & |
---|
1701 | indmin, indmax, & |
---|
1702 | pttab_child, petab_child, & |
---|
1703 | s_child, s_parent, & |
---|
1704 | ds_child, ds_parent ) |
---|
1705 | !--------------------------------------------------------------------------------------------------- |
---|
1706 | integer, dimension(5), intent(in) :: type_interp |
---|
1707 | integer, dimension(5), intent(in) :: indmin, indmax |
---|
1708 | integer, dimension(5), intent(in) :: pttab_child, petab_child |
---|
1709 | real, dimension(5), intent(in) :: s_child, s_parent |
---|
1710 | real, dimension(5), intent(in) :: ds_child, ds_parent |
---|
1711 | real, dimension( & |
---|
1712 | indmin(1):indmax(1), & |
---|
1713 | indmin(2):indmax(2), & |
---|
1714 | indmin(3):indmax(3), & |
---|
1715 | indmin(4):indmax(4), & |
---|
1716 | indmin(5):indmax(5)), intent(in) :: tabin |
---|
1717 | real, dimension( & |
---|
1718 | pttab_child(1):petab_child(1), & |
---|
1719 | pttab_child(2):petab_child(2), & |
---|
1720 | pttab_child(3):petab_child(3), & |
---|
1721 | pttab_child(4):petab_child(4), & |
---|
1722 | pttab_child(5):petab_child(5)), intent(out) :: tabout |
---|
1723 | !--------------------------------------------------------------------------------------------------- |
---|
1724 | real, dimension( & |
---|
1725 | pttab_child(1):petab_child(1), & |
---|
1726 | pttab_child(2):petab_child(2), & |
---|
1727 | pttab_child(3):petab_child(3), & |
---|
1728 | pttab_child(4):petab_child(4), & |
---|
1729 | indmin(5):indmax(5)) :: tabtemp |
---|
1730 | integer :: i, j, k, l, m |
---|
1731 | ! |
---|
1732 | do m = indmin(5), indmax(5) |
---|
1733 | call Agrif_Interp_4D_recursive(type_interp(1:4), & |
---|
1734 | tabin(indmin(1):indmax(1), & |
---|
1735 | indmin(2):indmax(2), & |
---|
1736 | indmin(3):indmax(3), & |
---|
1737 | indmin(4):indmax(4),m), & |
---|
1738 | tabtemp(pttab_child(1):petab_child(1), & |
---|
1739 | pttab_child(2):petab_child(2), & |
---|
1740 | pttab_child(3):petab_child(3), & |
---|
1741 | pttab_child(4):petab_child(4), m), & |
---|
1742 | indmin(1:4),indmax(1:4), & |
---|
1743 | pttab_child(1:4), petab_child(1:4), & |
---|
1744 | s_child(1:4), s_parent(1:4), & |
---|
1745 | ds_child(1:4), ds_parent(1:4) ) |
---|
1746 | enddo |
---|
1747 | ! |
---|
1748 | do l = pttab_child(4), petab_child(4) |
---|
1749 | do k = pttab_child(3), petab_child(3) |
---|
1750 | do j = pttab_child(2), petab_child(2) |
---|
1751 | do i = pttab_child(1), petab_child(1) |
---|
1752 | call Agrif_InterpBase(type_interp(5), & |
---|
1753 | tabtemp(i,j,k,l,indmin(5):indmax(5)), & |
---|
1754 | tabout(i,j,k,l,pttab_child(5):petab_child(5)), & |
---|
1755 | indmin(5), indmax(5), & |
---|
1756 | pttab_child(5), petab_child(5), & |
---|
1757 | s_parent(5), s_child(5), & |
---|
1758 | ds_parent(5), ds_child(5) ) |
---|
1759 | enddo |
---|
1760 | enddo |
---|
1761 | enddo |
---|
1762 | enddo |
---|
1763 | !--------------------------------------------------------------------------------------------------- |
---|
1764 | end subroutine Agrif_Interp_5D_recursive |
---|
1765 | !=================================================================================================== |
---|
1766 | ! |
---|
1767 | !=================================================================================================== |
---|
1768 | ! subroutine Agrif_Interp_6D_Recursive |
---|
1769 | ! |
---|
1770 | !> Subroutine for the interpolation of a 6D grid variable. |
---|
1771 | !> It calls #Agrif_Interp_5D_recursive and Agrif_InterpBase. |
---|
1772 | !--------------------------------------------------------------------------------------------------- |
---|
1773 | subroutine Agrif_Interp_6D_recursive ( type_interp, tabin, tabout, & |
---|
1774 | indmin, indmax, & |
---|
1775 | pttab_child, petab_child, & |
---|
1776 | s_child, s_parent, & |
---|
1777 | ds_child, ds_parent ) |
---|
1778 | !--------------------------------------------------------------------------------------------------- |
---|
1779 | integer, dimension(6), intent(in) :: type_interp |
---|
1780 | integer, dimension(6), intent(in) :: indmin, indmax |
---|
1781 | integer, dimension(6), intent(in) :: pttab_child, petab_child |
---|
1782 | real, dimension(6), intent(in) :: s_child, s_parent |
---|
1783 | real, dimension(6), intent(in) :: ds_child, ds_parent |
---|
1784 | real, dimension( & |
---|
1785 | indmin(1):indmax(1), & |
---|
1786 | indmin(2):indmax(2), & |
---|
1787 | indmin(3):indmax(3), & |
---|
1788 | indmin(4):indmax(4), & |
---|
1789 | indmin(5):indmax(5), & |
---|
1790 | indmin(6):indmax(6)), intent(in) :: tabin |
---|
1791 | real, dimension( & |
---|
1792 | pttab_child(1):petab_child(1), & |
---|
1793 | pttab_child(2):petab_child(2), & |
---|
1794 | pttab_child(3):petab_child(3), & |
---|
1795 | pttab_child(4):petab_child(4), & |
---|
1796 | pttab_child(5):petab_child(5), & |
---|
1797 | pttab_child(6):petab_child(6)), intent(out) :: tabout |
---|
1798 | !--------------------------------------------------------------------------------------------------- |
---|
1799 | real, dimension( & |
---|
1800 | pttab_child(1):petab_child(1), & |
---|
1801 | pttab_child(2):petab_child(2), & |
---|
1802 | pttab_child(3):petab_child(3), & |
---|
1803 | pttab_child(4):petab_child(4), & |
---|
1804 | pttab_child(5):petab_child(5), & |
---|
1805 | indmin(6):indmax(6)) :: tabtemp |
---|
1806 | integer :: i, j, k, l, m, n |
---|
1807 | ! |
---|
1808 | do n = indmin(6), indmax(6) |
---|
1809 | call Agrif_Interp_5D_recursive(type_interp(1:5), & |
---|
1810 | tabin(indmin(1):indmax(1), & |
---|
1811 | indmin(2):indmax(2), & |
---|
1812 | indmin(3):indmax(3), & |
---|
1813 | indmin(4):indmax(4), & |
---|
1814 | indmin(5):indmax(5), n), & |
---|
1815 | tabtemp(pttab_child(1):petab_child(1), & |
---|
1816 | pttab_child(2):petab_child(2), & |
---|
1817 | pttab_child(3):petab_child(3), & |
---|
1818 | pttab_child(4):petab_child(4), & |
---|
1819 | pttab_child(5):petab_child(5), n), & |
---|
1820 | indmin(1:5),indmax(1:5), & |
---|
1821 | pttab_child(1:5), petab_child(1:5), & |
---|
1822 | s_child(1:5), s_parent(1:5), & |
---|
1823 | ds_child(1:5),ds_parent(1:5) ) |
---|
1824 | enddo |
---|
1825 | ! |
---|
1826 | do m = pttab_child(5), petab_child(5) |
---|
1827 | do l = pttab_child(4), petab_child(4) |
---|
1828 | do k = pttab_child(3), petab_child(3) |
---|
1829 | do j = pttab_child(2), petab_child(2) |
---|
1830 | do i = pttab_child(1), petab_child(1) |
---|
1831 | call Agrif_InterpBase(type_interp(6), & |
---|
1832 | tabtemp(i,j,k,l,m,indmin(6):indmax(6)), & |
---|
1833 | tabout(i,j,k,l,m,pttab_child(6):petab_child(6)), & |
---|
1834 | indmin(6), indmax(6), & |
---|
1835 | pttab_child(6), petab_child(6), & |
---|
1836 | s_parent(6), s_child(6), & |
---|
1837 | ds_parent(6), ds_child(6) ) |
---|
1838 | enddo |
---|
1839 | enddo |
---|
1840 | enddo |
---|
1841 | enddo |
---|
1842 | enddo |
---|
1843 | !--------------------------------------------------------------------------------------------------- |
---|
1844 | end subroutine Agrif_Interp_6D_recursive |
---|
1845 | !=================================================================================================== |
---|
1846 | ! |
---|
1847 | !=================================================================================================== |
---|
1848 | ! subroutine Agrif_InterpBase |
---|
1849 | ! |
---|
1850 | !> Calls the interpolation method chosen by the user (linear, lagrange, spline, etc.). |
---|
1851 | !--------------------------------------------------------------------------------------------------- |
---|
1852 | subroutine Agrif_InterpBase ( type_interp, parenttab, childtab, indmin, indmax, & |
---|
1853 | pttab_child, petab_child, & |
---|
1854 | s_parent, s_child, ds_parent, ds_child ) |
---|
1855 | !--------------------------------------------------------------------------------------------------- |
---|
1856 | INTEGER :: type_interp |
---|
1857 | INTEGER :: indmin, indmax |
---|
1858 | INTEGER :: pttab_child, petab_child |
---|
1859 | REAL, DIMENSION(indmin:indmax), INTENT(IN) :: parenttab |
---|
1860 | REAL, DIMENSION(pttab_child:petab_child), INTENT(OUT) :: childtab |
---|
1861 | REAL :: s_parent, s_child |
---|
1862 | REAL :: ds_parent,ds_child |
---|
1863 | ! |
---|
1864 | if ( (indmin == indmax) .and. (pttab_child == petab_child) ) then |
---|
1865 | ! |
---|
1866 | childtab(pttab_child) = parenttab(indmin) |
---|
1867 | ! |
---|
1868 | elseif (type_interp == Agrif_LINEAR) then ! Linear interpolation |
---|
1869 | ! |
---|
1870 | call Agrif_basicinterp_linear1D(parenttab,childtab, & |
---|
1871 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1872 | s_parent,s_child,ds_parent,ds_child) |
---|
1873 | ! |
---|
1874 | elseif ( type_interp == Agrif_PPM ) then ! PPM interpolation |
---|
1875 | |
---|
1876 | call PPM1d(parenttab,childtab, & |
---|
1877 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1878 | s_parent,s_child,ds_parent,ds_child) |
---|
1879 | ! |
---|
1880 | elseif ( type_interp == Agrif_PPM_LIM ) then ! PPM interpolation |
---|
1881 | |
---|
1882 | call PPM1d_lim(parenttab,childtab, & |
---|
1883 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1884 | s_parent,s_child,ds_parent,ds_child) |
---|
1885 | ! |
---|
1886 | elseif (type_interp == Agrif_LAGRANGE) then ! Lagrange interpolation |
---|
1887 | ! |
---|
1888 | call lagrange1D(parenttab,childtab, & |
---|
1889 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1890 | s_parent,s_child,ds_parent,ds_child) |
---|
1891 | ! |
---|
1892 | elseif (type_interp == Agrif_ENO) then ! Eno interpolation |
---|
1893 | ! |
---|
1894 | call ENO1d(parenttab,childtab, & |
---|
1895 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1896 | s_parent,s_child,ds_parent,ds_child) |
---|
1897 | ! |
---|
1898 | elseif (type_interp == Agrif_WENO) then ! Weno interpolation |
---|
1899 | ! |
---|
1900 | call WENO1d(parenttab,childtab, & |
---|
1901 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1902 | s_parent,s_child,ds_parent,ds_child) |
---|
1903 | ! |
---|
1904 | elseif (type_interp == Agrif_LINEARCONSERV) then ! Linear conservative interpolation |
---|
1905 | ! |
---|
1906 | call Linear1dConserv(parenttab,childtab, & |
---|
1907 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1908 | s_parent,s_child,ds_parent,ds_child) |
---|
1909 | ! |
---|
1910 | elseif (type_interp == Agrif_LINEARCONSERVLIM) then !Linear conservative interpolation |
---|
1911 | ! |
---|
1912 | call Linear1dConservLim(parenttab,childtab, & |
---|
1913 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1914 | s_parent,s_child,ds_parent,ds_child) |
---|
1915 | ! |
---|
1916 | elseif (type_interp == Agrif_CONSTANT) then |
---|
1917 | ! |
---|
1918 | call Constant1d(parenttab,childtab, & |
---|
1919 | indmax-indmin+1,petab_child-pttab_child+1, & |
---|
1920 | s_parent,s_child,ds_parent,ds_child) |
---|
1921 | ! |
---|
1922 | endif |
---|
1923 | !--------------------------------------------------------------------------------------------------- |
---|
1924 | end subroutine Agrif_InterpBase |
---|
1925 | !=================================================================================================== |
---|
1926 | ! |
---|
1927 | !=================================================================================================== |
---|
1928 | ! subroutine Agrif_Find_list_interp |
---|
1929 | !--------------------------------------------------------------------------------------------------- |
---|
1930 | function Agrif_Find_list_interp ( list_interp, pttab, petab, pttab_Child, pttab_Parent, & |
---|
1931 | nbdim, indmin, indmax, indminglob, indmaxglob, & |
---|
1932 | pttruetab, cetruetab, memberin, & |
---|
1933 | parentarray_chunk,parentarray_chunk_decal,decal_chunks, & |
---|
1934 | correction_required,member_chuncks,nb_chunks & |
---|
1935 | #if defined AGRIF_MPI |
---|
1936 | ,indminglob2, indmaxglob2, parentarray, & |
---|
1937 | member, tab4t, memberinall, sendtoproc1, recvfromproc1 & |
---|
1938 | #endif |
---|
1939 | ) result(find_list_interp) |
---|
1940 | !--------------------------------------------------------------------------------------------------- |
---|
1941 | type(Agrif_List_Interp_Loc), pointer :: list_interp |
---|
1942 | integer, intent(in) :: nbdim |
---|
1943 | integer, dimension(nbdim), intent(in) :: pttab, petab, pttab_Child, pttab_Parent |
---|
1944 | integer, dimension(nbdim), intent(out) :: indmin, indmax |
---|
1945 | integer, dimension(nbdim), intent(out) :: indminglob, indmaxglob |
---|
1946 | integer, dimension(nbdim), intent(out) :: pttruetab, cetruetab |
---|
1947 | logical, intent(out) :: memberin |
---|
1948 | integer :: nb_chunks |
---|
1949 | integer, dimension(:,:,:,:), allocatable :: parentarray_chunk |
---|
1950 | integer, dimension(:,:,:,:), allocatable :: parentarray_chunk_decal |
---|
1951 | integer, dimension(:,:),allocatable :: decal_chunks |
---|
1952 | logical, dimension(:),allocatable :: correction_required |
---|
1953 | logical, dimension(:),allocatable :: member_chuncks |
---|
1954 | #if defined AGRIF_MPI |
---|
1955 | integer, dimension(nbdim), intent(out) :: indminglob2, indmaxglob2 |
---|
1956 | integer, dimension(nbdim,2,2), intent(out) :: parentarray |
---|
1957 | logical, intent(out) :: member |
---|
1958 | integer, dimension(nbdim,0:Agrif_Nbprocs-1,8), intent(out) :: tab4t |
---|
1959 | logical, dimension(0:Agrif_Nbprocs-1), intent(out) :: memberinall |
---|
1960 | logical, dimension(0:Agrif_Nbprocs-1), intent(out) :: sendtoproc1, recvfromproc1 |
---|
1961 | #endif |
---|
1962 | logical :: find_list_interp |
---|
1963 | ! |
---|
1964 | integer :: i |
---|
1965 | type(Agrif_List_Interp_Loc), pointer :: parcours |
---|
1966 | type(Agrif_Interp_Loc), pointer :: pil |
---|
1967 | |
---|
1968 | find_list_interp = .false. |
---|
1969 | |
---|
1970 | if ( .not. associated(list_interp) ) return |
---|
1971 | |
---|
1972 | parcours => list_interp |
---|
1973 | find_loop : do while ( associated(parcours) ) |
---|
1974 | |
---|
1975 | pil => parcours % interp_loc |
---|
1976 | |
---|
1977 | do i = 1,nbdim |
---|
1978 | if ( (pttab(i) /= pil % pttab(i)) .or. & |
---|
1979 | (petab(i) /= pil % petab(i)) .or. & |
---|
1980 | (pttab_child(i) /= pil % pttab_child(i)) .or. & |
---|
1981 | (pttab_parent(i) /= pil % pttab_parent(i)) ) then |
---|
1982 | parcours => parcours % suiv |
---|
1983 | cycle find_loop |
---|
1984 | endif |
---|
1985 | enddo |
---|
1986 | |
---|
1987 | indmin = pil % indmin(1:nbdim) |
---|
1988 | indmax = pil % indmax(1:nbdim) |
---|
1989 | |
---|
1990 | pttruetab = pil % pttruetab(1:nbdim) |
---|
1991 | cetruetab = pil % cetruetab(1:nbdim) |
---|
1992 | |
---|
1993 | #if !defined AGRIF_MPI |
---|
1994 | indminglob = pil % indminglob(1:nbdim) |
---|
1995 | indmaxglob = pil % indmaxglob(1:nbdim) |
---|
1996 | #else |
---|
1997 | indminglob = pil % indminglob2(1:nbdim) |
---|
1998 | indmaxglob = pil % indmaxglob2(1:nbdim) |
---|
1999 | indminglob2 = pil % indminglob2(1:nbdim) |
---|
2000 | indmaxglob2 = pil % indmaxglob2(1:nbdim) |
---|
2001 | parentarray = pil % parentarray(1:nbdim,:,:) |
---|
2002 | member = pil % member |
---|
2003 | tab4t = pil % tab4t(1:nbdim, 0:Agrif_Nbprocs-1, 1:8) |
---|
2004 | memberinall = pil % memberinall(0:Agrif_Nbprocs-1) |
---|
2005 | sendtoproc1 = pil % sendtoproc1(0:Agrif_Nbprocs-1) |
---|
2006 | recvfromproc1 = pil % recvfromproc1(0:Agrif_Nbprocs-1) |
---|
2007 | #endif |
---|
2008 | memberin = pil % memberin |
---|
2009 | |
---|
2010 | ! chunks |
---|
2011 | nb_chunks = pil % nb_chunks |
---|
2012 | Allocate(parentarray_chunk(nb_chunks,nbdim,2,2)) |
---|
2013 | parentarray_chunk = pil % parentarray_chunk |
---|
2014 | Allocate(parentarray_chunk_decal(nb_chunks,nbdim,2,2)) |
---|
2015 | parentarray_chunk_decal = pil % parentarray_chunk_decal |
---|
2016 | Allocate(correction_required(nb_chunks)) |
---|
2017 | correction_required = pil % correction_required |
---|
2018 | Allocate(decal_chunks(nb_chunks,nbdim)) |
---|
2019 | decal_chunks = pil % decal_chunks |
---|
2020 | Allocate(member_chuncks(nb_chunks)) |
---|
2021 | member_chuncks = pil % member_chuncks |
---|
2022 | |
---|
2023 | find_list_interp = .true. |
---|
2024 | exit find_loop |
---|
2025 | enddo find_loop |
---|
2026 | !--------------------------------------------------------------------------------------------------- |
---|
2027 | end function Agrif_Find_list_interp |
---|
2028 | !=================================================================================================== |
---|
2029 | ! |
---|
2030 | !=================================================================================================== |
---|
2031 | ! subroutine Agrif_AddTo_list_interp |
---|
2032 | !--------------------------------------------------------------------------------------------------- |
---|
2033 | subroutine Agrif_AddTo_list_interp ( list_interp, pttab, petab, pttab_Child, pttab_Parent, & |
---|
2034 | indmin, indmax, indminglob, indmaxglob, & |
---|
2035 | pttruetab, cetruetab, & |
---|
2036 | memberin, nbdim, & |
---|
2037 | parentarray_chunk,parentarray_chunk_decal,decal_chunks, & |
---|
2038 | correction_required,member_chuncks,nb_chunks & |
---|
2039 | #if defined AGRIF_MPI |
---|
2040 | ,indminglob2, indmaxglob2, & |
---|
2041 | parentarray, & |
---|
2042 | member, & |
---|
2043 | tab4t, memberinall, sendtoproc1, recvfromproc1 & |
---|
2044 | #endif |
---|
2045 | ) |
---|
2046 | !--------------------------------------------------------------------------------------------------- |
---|
2047 | type(Agrif_List_Interp_Loc), pointer :: list_interp |
---|
2048 | integer :: nbdim |
---|
2049 | integer, dimension(nbdim) :: pttab, petab, pttab_Child, pttab_Parent |
---|
2050 | integer, dimension(nbdim) :: indmin,indmax |
---|
2051 | integer, dimension(nbdim) :: indminglob, indmaxglob |
---|
2052 | integer, dimension(nbdim) :: pttruetab, cetruetab |
---|
2053 | logical :: memberin |
---|
2054 | integer :: nb_chunks |
---|
2055 | integer, dimension(:,:,:,:), allocatable :: parentarray_chunk |
---|
2056 | integer, dimension(:,:,:,:), allocatable :: parentarray_chunk_decal |
---|
2057 | integer, dimension(:,:),allocatable :: decal_chunks |
---|
2058 | logical, dimension(:),allocatable :: correction_required |
---|
2059 | logical, dimension(:),allocatable :: member_chuncks |
---|
2060 | #if defined AGRIF_MPI |
---|
2061 | integer, dimension(nbdim,2,2) :: parentarray |
---|
2062 | logical :: member |
---|
2063 | integer, dimension(nbdim) :: indminglob2,indmaxglob2 |
---|
2064 | integer, dimension(nbdim,0:Agrif_Nbprocs-1,8) :: tab4t |
---|
2065 | logical, dimension(0:Agrif_Nbprocs-1) :: memberinall |
---|
2066 | logical, dimension(0:Agrif_Nbprocs-1) :: sendtoproc1 |
---|
2067 | logical, dimension(0:Agrif_Nbprocs-1) :: recvfromproc1 |
---|
2068 | #endif |
---|
2069 | ! |
---|
2070 | type(Agrif_List_Interp_Loc), pointer :: parcours |
---|
2071 | type(Agrif_Interp_Loc), pointer :: pil |
---|
2072 | ! |
---|
2073 | allocate(parcours) |
---|
2074 | allocate(parcours % interp_loc) |
---|
2075 | |
---|
2076 | pil => parcours % interp_loc |
---|
2077 | |
---|
2078 | pil % pttab(1:nbdim) = pttab(1:nbdim) |
---|
2079 | pil % petab(1:nbdim) = petab(1:nbdim) |
---|
2080 | pil % pttab_child(1:nbdim) = pttab_child(1:nbdim) |
---|
2081 | pil % pttab_parent(1:nbdim) = pttab_parent(1:nbdim) |
---|
2082 | |
---|
2083 | pil % indmin(1:nbdim) = indmin(1:nbdim) |
---|
2084 | pil % indmax(1:nbdim) = indmax(1:nbdim) |
---|
2085 | |
---|
2086 | pil % memberin = memberin |
---|
2087 | #if !defined AGRIF_MPI |
---|
2088 | pil % indminglob(1:nbdim) = indminglob(1:nbdim) |
---|
2089 | pil % indmaxglob(1:nbdim) = indmaxglob(1:nbdim) |
---|
2090 | #else |
---|
2091 | pil % indminglob2(1:nbdim) = indminglob2(1:nbdim) |
---|
2092 | pil % indmaxglob2(1:nbdim) = indmaxglob2(1:nbdim) |
---|
2093 | pil % parentarray(1:nbdim,:,:) = parentarray(1:nbdim,:,:) |
---|
2094 | pil % member = member |
---|
2095 | allocate(pil % tab4t(nbdim, 0:Agrif_Nbprocs-1, 8)) |
---|
2096 | allocate(pil % memberinall(0:Agrif_Nbprocs-1)) |
---|
2097 | allocate(pil % sendtoproc1(0:Agrif_Nbprocs-1)) |
---|
2098 | allocate(pil % recvfromproc1(0:Agrif_Nbprocs-1)) |
---|
2099 | pil % tab4t = tab4t |
---|
2100 | pil % memberinall = memberinall |
---|
2101 | pil % sendtoproc1 = sendtoproc1 |
---|
2102 | pil % recvfromproc1 = recvfromproc1 |
---|
2103 | #endif |
---|
2104 | |
---|
2105 | pil % pttruetab(1:nbdim) = pttruetab(1:nbdim) |
---|
2106 | pil % cetruetab(1:nbdim) = cetruetab(1:nbdim) |
---|
2107 | |
---|
2108 | ! chunks |
---|
2109 | pil % nb_chunks = nb_chunks |
---|
2110 | allocate(pil % parentarray_chunk(nb_chunks,nbdim,2,2)) |
---|
2111 | allocate(pil % parentarray_chunk_decal(nb_chunks,nbdim,2,2)) |
---|
2112 | allocate(pil % correction_required(nb_chunks)) |
---|
2113 | allocate(pil % decal_chunks(nb_chunks,nbdim)) |
---|
2114 | allocate(pil % member_chuncks(nb_chunks)) |
---|
2115 | |
---|
2116 | pil % parentarray_chunk = parentarray_chunk |
---|
2117 | pil % parentarray_chunk_decal = parentarray_chunk_decal |
---|
2118 | pil % correction_required = correction_required |
---|
2119 | pil % decal_chunks = decal_chunks |
---|
2120 | pil % member_chuncks = member_chuncks |
---|
2121 | |
---|
2122 | |
---|
2123 | parcours % suiv => list_interp |
---|
2124 | list_interp => parcours |
---|
2125 | !--------------------------------------------------------------------------------------------------- |
---|
2126 | end subroutine Agrif_Addto_list_interp |
---|
2127 | !=================================================================================================== |
---|
2128 | ! |
---|
2129 | end module Agrif_Interpolation |
---|