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