1 | ! |
---|
2 | ! $Id: modupdate.F 779 2007-12-22 17:04:17Z rblod $ |
---|
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_Update |
---|
26 | C |
---|
27 | Module Agrif_Update |
---|
28 | C |
---|
29 | CCC Description: |
---|
30 | CCC Module to update a parent grid from its child grids |
---|
31 | C |
---|
32 | C Modules used: |
---|
33 | C |
---|
34 | Use Agrif_Updatebasic |
---|
35 | c Use Agrif_Boundary |
---|
36 | Use Agrif_Arrays |
---|
37 | Use Agrif_CurgridFunctions |
---|
38 | Use Agrif_Mask |
---|
39 | #ifdef 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_Update_1d |
---|
53 | C ************************************************************************** |
---|
54 | C |
---|
55 | Subroutine Agrif_Update_1d(TypeUpdate,parent,child,tab,deb,fin, |
---|
56 | & procname) |
---|
57 | C |
---|
58 | CCC Description: |
---|
59 | CCC Subroutine to update a 1D grid variable on the parent grid. |
---|
60 | C |
---|
61 | C Declarations: |
---|
62 | C |
---|
63 | |
---|
64 | C |
---|
65 | C Arguments |
---|
66 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update (copy or average) |
---|
67 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
68 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
69 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
70 | INTEGER, DIMENSION(6) :: deb,fin ! Positions where interpolations |
---|
71 | ! are done on the fine grid |
---|
72 | External :: procname |
---|
73 | Optional :: procname |
---|
74 | REAL, DIMENSION( |
---|
75 | & child%var%lb(1):child%var%ub(1) |
---|
76 | & ), Target :: tab ! Result |
---|
77 | C |
---|
78 | C |
---|
79 | C Definition of a temporary AGRIF_PVariable data TYPE |
---|
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 = 1 |
---|
87 | C |
---|
88 | C Values on the current grid used for the update |
---|
89 | C childtemp % var % array1 => tab |
---|
90 | |
---|
91 | childtemp % var % lb = child % var % lb |
---|
92 | childtemp % var % ub = child % var % ub |
---|
93 | |
---|
94 | C childtemp % var % list_update => child%var%list_update |
---|
95 | |
---|
96 | C |
---|
97 | |
---|
98 | IF (present(procname)) THEN |
---|
99 | CALL Agrif_UpdateVariable |
---|
100 | & (TypeUpdate,parent,child,deb,fin,procname) |
---|
101 | ELSE |
---|
102 | CALL Agrif_UpdateVariable |
---|
103 | & (TypeUpdate,parent,child,deb,fin) |
---|
104 | ENDIF |
---|
105 | C |
---|
106 | C child % var % list_update => childtemp%var%list_update |
---|
107 | |
---|
108 | deallocate(childtemp % var) |
---|
109 | C |
---|
110 | C |
---|
111 | End Subroutine Agrif_Update_1D |
---|
112 | C |
---|
113 | C |
---|
114 | C |
---|
115 | C ************************************************************************** |
---|
116 | CCC Subroutine Agrif_Update_2d |
---|
117 | C ************************************************************************** |
---|
118 | C |
---|
119 | |
---|
120 | Subroutine Agrif_Update_2d(TypeUpdate,parent,child,tab,deb,fin, |
---|
121 | & procname) |
---|
122 | C |
---|
123 | CCC Description: |
---|
124 | CCC Subroutine to update a 2D grid variable on the parent grid. |
---|
125 | C |
---|
126 | C Declarations: |
---|
127 | C |
---|
128 | |
---|
129 | C |
---|
130 | C Arguments |
---|
131 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update (copy or average) |
---|
132 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
133 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
134 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
135 | INTEGER, DIMENSION(6) :: deb,fin ! Positions where interpolations |
---|
136 | ! are done on the fine grid |
---|
137 | |
---|
138 | External :: procname |
---|
139 | Optional :: procname |
---|
140 | |
---|
141 | REAL, DIMENSION( |
---|
142 | & child%var%lb(1):child%var%ub(1), |
---|
143 | & child%var%lb(2):child%var%ub(2) |
---|
144 | & ), Target :: tab ! Result |
---|
145 | C |
---|
146 | C |
---|
147 | C Definition of a temporary AGRIF_PVariable data TYPE |
---|
148 | allocate(childtemp % var) |
---|
149 | C |
---|
150 | C Pointer on the root variable |
---|
151 | childtemp % var % root_var => child % var %root_var |
---|
152 | C |
---|
153 | C Number of dimensions of the grid variable |
---|
154 | childtemp % var % nbdim = 2 |
---|
155 | C |
---|
156 | C Values on the current grid used for the update |
---|
157 | C childtemp % var % array2 => tab |
---|
158 | |
---|
159 | childtemp % var % lb = child % var % lb |
---|
160 | childtemp % var % ub = child % var % ub |
---|
161 | |
---|
162 | C childtemp % var % list_update => child%var%list_update |
---|
163 | C |
---|
164 | IF (present(procname)) THEN |
---|
165 | CALL Agrif_UpdateVariable |
---|
166 | & (TypeUpdate,parent,child,deb,fin,procname) |
---|
167 | ELSE |
---|
168 | CALL Agrif_UpdateVariable |
---|
169 | & (TypeUpdate,parent,child,deb,fin) |
---|
170 | ENDIF |
---|
171 | C |
---|
172 | C child % var % list_update => childtemp%var%list_update |
---|
173 | |
---|
174 | deallocate(childtemp % var) |
---|
175 | C |
---|
176 | C |
---|
177 | End Subroutine Agrif_Update_2D |
---|
178 | C |
---|
179 | C |
---|
180 | C |
---|
181 | C ************************************************************************** |
---|
182 | CCC Subroutine Agrif_Update_3d |
---|
183 | C ************************************************************************** |
---|
184 | C |
---|
185 | Subroutine Agrif_Update_3d(TypeUpdate,parent,child,tab,deb,fin, |
---|
186 | & procname) |
---|
187 | C |
---|
188 | CCC Description: |
---|
189 | CCC Subroutine to update a 3D grid variable on the parent grid. |
---|
190 | C |
---|
191 | C Declarations: |
---|
192 | C |
---|
193 | |
---|
194 | C |
---|
195 | C Arguments |
---|
196 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update (copy or average) |
---|
197 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
198 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
199 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
200 | INTEGER, DIMENSION(6) :: deb,fin ! Positions where interpolations |
---|
201 | ! are done on the fine grid |
---|
202 | External :: procname |
---|
203 | Optional :: procname |
---|
204 | |
---|
205 | REAL, DIMENSION( |
---|
206 | & child%var%lb(1):child%var%ub(1), |
---|
207 | & child%var%lb(2):child%var%ub(2), |
---|
208 | & child%var%lb(3):child%var%ub(3) |
---|
209 | & ), Target :: tab ! Results |
---|
210 | C |
---|
211 | C |
---|
212 | C Definition of a temporary AGRIF_PVariable data TYPE |
---|
213 | allocate(childtemp % var) |
---|
214 | C |
---|
215 | C Pointer on the root variable |
---|
216 | childtemp % var % root_var => child % var %root_var |
---|
217 | C |
---|
218 | C Number of dimensions of the grid variable |
---|
219 | childtemp % var % nbdim = 3 |
---|
220 | C |
---|
221 | C Values on the current grid used for the update |
---|
222 | C childtemp % var % array3 => tab |
---|
223 | |
---|
224 | childtemp % var % lb = child % var % lb |
---|
225 | childtemp % var % ub = child % var % ub |
---|
226 | |
---|
227 | |
---|
228 | C childtemp % var % list_update => child%var%list_update |
---|
229 | C |
---|
230 | IF (present(procname)) THEN |
---|
231 | CALL Agrif_UpdateVariable |
---|
232 | & (TypeUpdate,parent,child,deb,fin,procname) |
---|
233 | ELSE |
---|
234 | CALL Agrif_UpdateVariable |
---|
235 | & (TypeUpdate,parent,child,deb,fin) |
---|
236 | ENDIF |
---|
237 | C |
---|
238 | C child % var % list_update => childtemp%var%list_update |
---|
239 | |
---|
240 | DEALLOCATE(childtemp % var) |
---|
241 | C |
---|
242 | C |
---|
243 | End Subroutine Agrif_Update_3D |
---|
244 | C |
---|
245 | C |
---|
246 | C |
---|
247 | C ************************************************************************** |
---|
248 | CCC Subroutine Agrif_Update_4d |
---|
249 | C ************************************************************************** |
---|
250 | C |
---|
251 | Subroutine Agrif_Update_4d(TypeUpdate,parent,child,tab,deb,fin, |
---|
252 | & procname) |
---|
253 | C |
---|
254 | CCC Description: |
---|
255 | CCC Subroutine to update a 4D grid variable on the parent grid. |
---|
256 | C |
---|
257 | C Declarations: |
---|
258 | C |
---|
259 | |
---|
260 | C |
---|
261 | C Arguments |
---|
262 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update (copy or average) |
---|
263 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
264 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
265 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
266 | INTEGER, DIMENSION(6) :: deb,fin ! Positions where interpolations |
---|
267 | ! are done on the fine grid |
---|
268 | External :: procname |
---|
269 | Optional :: procname |
---|
270 | REAL, DIMENSION( |
---|
271 | & child%var%lb(1):child%var%ub(1), |
---|
272 | & child%var%lb(2):child%var%ub(2), |
---|
273 | & child%var%lb(3):child%var%ub(3), |
---|
274 | & child%var%lb(4):child%var%ub(4) |
---|
275 | & ), Target :: tab ! Results |
---|
276 | C |
---|
277 | C |
---|
278 | C Definition of a temporary AGRIF_PVariable data TYPE |
---|
279 | allocate(childtemp % var) |
---|
280 | C |
---|
281 | C Pointer on the root variable |
---|
282 | childtemp % var % root_var => child % var %root_var |
---|
283 | C |
---|
284 | C Number of dimensions of the grid variable |
---|
285 | childtemp % var % nbdim = 4 |
---|
286 | C |
---|
287 | C Values on the current grid used for the update |
---|
288 | C childtemp % var % array4 => tab |
---|
289 | |
---|
290 | childtemp % var % lb = child % var % lb |
---|
291 | childtemp % var % ub = child % var % ub |
---|
292 | |
---|
293 | |
---|
294 | C childtemp % var % list_update => child%var%list_update |
---|
295 | |
---|
296 | C |
---|
297 | IF (present(procname)) THEN |
---|
298 | CALL Agrif_UpdateVariable |
---|
299 | & (TypeUpdate,parent,child,deb,fin,procname) |
---|
300 | ELSE |
---|
301 | CALL Agrif_UpdateVariable |
---|
302 | & (TypeUpdate,parent,child,deb,fin) |
---|
303 | ENDIF |
---|
304 | |
---|
305 | C child % var % list_update => childtemp%var%list_update |
---|
306 | C |
---|
307 | deallocate(childtemp % var) |
---|
308 | C |
---|
309 | C |
---|
310 | End Subroutine Agrif_Update_4D |
---|
311 | C |
---|
312 | C |
---|
313 | C |
---|
314 | C ************************************************************************** |
---|
315 | CCC Subroutine Agrif_Update_5d |
---|
316 | C ************************************************************************** |
---|
317 | C |
---|
318 | Subroutine Agrif_Update_5d(TypeUpdate,parent,child,tab,deb,fin, |
---|
319 | & procname) |
---|
320 | C |
---|
321 | CCC Description: |
---|
322 | CCC Subroutine to update a 5D grid variable on the parent grid. |
---|
323 | C |
---|
324 | C Declarations: |
---|
325 | C |
---|
326 | |
---|
327 | C |
---|
328 | C Arguments |
---|
329 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update (copy or average) |
---|
330 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
331 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
332 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
333 | INTEGER, DIMENSION(6) :: deb,fin ! Positions where interpolations |
---|
334 | ! are done on the fine grid |
---|
335 | External :: procname |
---|
336 | Optional :: procname |
---|
337 | |
---|
338 | REAL, DIMENSION( |
---|
339 | & child%var%lb(1):child%var%ub(1), |
---|
340 | & child%var%lb(2):child%var%ub(2), |
---|
341 | & child%var%lb(3):child%var%ub(3), |
---|
342 | & child%var%lb(4):child%var%ub(4), |
---|
343 | & child%var%lb(5):child%var%ub(5) |
---|
344 | & ), Target :: tab ! Results |
---|
345 | C |
---|
346 | C |
---|
347 | C Definition of a temporary AGRIF_PVariable data TYPE |
---|
348 | allocate(childtemp % var) |
---|
349 | C |
---|
350 | C Pointer on the root variable |
---|
351 | childtemp % var % root_var => child % var %root_var |
---|
352 | C |
---|
353 | C Number of dimensions of the grid variable |
---|
354 | childtemp % var % nbdim = 5 |
---|
355 | C |
---|
356 | C Values on the current grid used for the update |
---|
357 | C childtemp % var % array5 => tab |
---|
358 | |
---|
359 | childtemp % var % lb = child % var % lb |
---|
360 | childtemp % var % ub = child % var % ub |
---|
361 | |
---|
362 | C childtemp % var % list_update => child%var%list_update |
---|
363 | C |
---|
364 | IF (present(procname)) THEN |
---|
365 | CALL Agrif_UpdateVariable |
---|
366 | & (TypeUpdate,parent,child,deb,fin,procname) |
---|
367 | ELSE |
---|
368 | CALL Agrif_UpdateVariable |
---|
369 | & (TypeUpdate,parent,child,deb,fin) |
---|
370 | ENDIF |
---|
371 | |
---|
372 | C child % var % list_update => childtemp%var%list_update |
---|
373 | |
---|
374 | C |
---|
375 | deallocate(childtemp % var) |
---|
376 | C |
---|
377 | C |
---|
378 | End Subroutine Agrif_Update_5D |
---|
379 | C |
---|
380 | C |
---|
381 | C |
---|
382 | C |
---|
383 | C ************************************************************************** |
---|
384 | CCC Subroutine Agrif_Update_6d |
---|
385 | C ************************************************************************** |
---|
386 | C |
---|
387 | Subroutine Agrif_Update_6d(TypeUpdate,parent,child,tab,deb,fin) |
---|
388 | C |
---|
389 | CCC Description: |
---|
390 | CCC Subroutine to update a 6D grid variable on the parent grid. |
---|
391 | C |
---|
392 | C Declarations: |
---|
393 | C |
---|
394 | |
---|
395 | C |
---|
396 | C Arguments |
---|
397 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update (copy or average) |
---|
398 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
399 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
400 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
401 | INTEGER, DIMENSION(6) :: deb,fin ! Positions where interpolations |
---|
402 | ! are done on the fine grid |
---|
403 | REAL, DIMENSION( |
---|
404 | & child%var%lb(1):child%var%ub(1), |
---|
405 | & child%var%lb(2):child%var%ub(2), |
---|
406 | & child%var%lb(3):child%var%ub(3), |
---|
407 | & child%var%lb(4):child%var%ub(4), |
---|
408 | & child%var%lb(5):child%var%ub(5), |
---|
409 | & child%var%lb(6):child%var%ub(6) |
---|
410 | & ), Target :: tab ! Results |
---|
411 | C |
---|
412 | C |
---|
413 | C Definition of a temporary AGRIF_PVariable data TYPE |
---|
414 | allocate(childtemp % var) |
---|
415 | C |
---|
416 | C Pointer on the root variable |
---|
417 | childtemp % var % root_var => child % var %root_var |
---|
418 | C |
---|
419 | C Number of dimensions of the grid variable |
---|
420 | childtemp % var % nbdim = 6 |
---|
421 | C |
---|
422 | C Values on the current grid used for the update |
---|
423 | C childtemp % var % array6 => tab |
---|
424 | |
---|
425 | childtemp % var % lb = child % var % lb |
---|
426 | childtemp % var % ub = child % var % ub |
---|
427 | |
---|
428 | C childtemp % var % list_update => child%var%list_update |
---|
429 | C |
---|
430 | Call Agrif_UpdateVariable |
---|
431 | & (TypeUpdate,parent,child,deb,fin) |
---|
432 | |
---|
433 | C child % var % list_update => childtemp%var%list_update |
---|
434 | |
---|
435 | C |
---|
436 | deallocate(childtemp % var) |
---|
437 | C |
---|
438 | C |
---|
439 | End Subroutine Agrif_Update_6D |
---|
440 | C |
---|
441 | C |
---|
442 | C |
---|
443 | C ************************************************************************** |
---|
444 | C Subroutine Agrif_UpdateVariable |
---|
445 | C ************************************************************************** |
---|
446 | C |
---|
447 | Subroutine Agrif_UpdateVariable(TypeUpdate,parent,child,deb,fin, |
---|
448 | & procname) |
---|
449 | C |
---|
450 | CCC Description: |
---|
451 | CCC Subroutine to set arguments of Agrif_UpdatenD, n being the number of |
---|
452 | C dimensions of the grid variable. |
---|
453 | C |
---|
454 | CC Declarations: |
---|
455 | C |
---|
456 | c |
---|
457 | C |
---|
458 | C Scalar argument |
---|
459 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update (copy or average) |
---|
460 | C Data TYPE arguments |
---|
461 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
462 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
463 | INTEGER, DIMENSION(6) :: deb,fin ! Positions where interpolations |
---|
464 | ! are calculated |
---|
465 | External :: procname |
---|
466 | Optional :: procname |
---|
467 | |
---|
468 | C |
---|
469 | C Local scalars |
---|
470 | INTEGER :: nbdim ! Number of dimensions of the current |
---|
471 | ! grid |
---|
472 | INTEGER ,DIMENSION(6) :: pttab_child |
---|
473 | INTEGER ,DIMENSION(6) :: petab_child |
---|
474 | INTEGER ,DIMENSION(6) :: pttab_parent |
---|
475 | REAL ,DIMENSION(6) :: s_child,s_parent |
---|
476 | REAL ,DIMENSION(6) :: ds_child,ds_parent |
---|
477 | INTEGER,DIMENSION(6) :: loctab_Child ! Indicates if the child |
---|
478 | ! grid has a common border with |
---|
479 | ! the root grid |
---|
480 | TYPE(AGRIF_Variable), Pointer :: root ! Variable on the root grid |
---|
481 | INTEGER,DIMENSION(6) :: posvartab_Child ! Position of the |
---|
482 | ! variable on the cell |
---|
483 | INTEGER,DIMENSION(6) :: nbtab_Child ! Number of the cells |
---|
484 | INTEGER :: n |
---|
485 | LOGICAL :: wholeupdate |
---|
486 | C |
---|
487 | C |
---|
488 | |
---|
489 | loctab_child(:) = 0 |
---|
490 | C |
---|
491 | root => child % var % root_var |
---|
492 | nbdim = root % nbdim |
---|
493 | C |
---|
494 | do n = 1,nbdim |
---|
495 | posvartab_child(n) = root % posvar(n) |
---|
496 | enddo |
---|
497 | C |
---|
498 | |
---|
499 | Call PreProcessToInterpOrUpdate(parent,child, |
---|
500 | & petab_Child(1:nbdim), |
---|
501 | & pttab_Child(1:nbdim),pttab_Parent(1:nbdim), |
---|
502 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
503 | & ds_Child(1:nbdim),ds_Parent(1:nbdim), |
---|
504 | & nbdim) |
---|
505 | C |
---|
506 | C |
---|
507 | do n = 1,nbdim |
---|
508 | C |
---|
509 | Select case(root % interptab(n)) |
---|
510 | C |
---|
511 | case('x') ! x DIMENSION |
---|
512 | C |
---|
513 | nbtab_Child(n) = Agrif_Curgrid % nb(1) |
---|
514 | C |
---|
515 | case('y') ! y DIMENSION |
---|
516 | C |
---|
517 | nbtab_Child(n) = Agrif_Curgrid % nb(2) |
---|
518 | C |
---|
519 | case('z') ! z DIMENSION |
---|
520 | C |
---|
521 | nbtab_Child(n) = Agrif_Curgrid % nb(3) |
---|
522 | C |
---|
523 | case('N') ! No space DIMENSION |
---|
524 | C |
---|
525 | |
---|
526 | nbtab_Child(n) = child % var % ub(n) - child % var % lb(n) |
---|
527 | C |
---|
528 | C No interpolation but only a copy of the values of the grid variable |
---|
529 | C |
---|
530 | posvartab_child(n) = 1 |
---|
531 | |
---|
532 | loctab_child(n) = -3 |
---|
533 | C |
---|
534 | End select |
---|
535 | C |
---|
536 | enddo |
---|
537 | |
---|
538 | C Call to a procedure of update according to the number of dimensions of |
---|
539 | C the grid variable |
---|
540 | |
---|
541 | wholeupdate = .FALSE. |
---|
542 | |
---|
543 | do n=1,nbdim |
---|
544 | if (loctab_child(n) /= -3) then |
---|
545 | if (deb(n)>fin(n)) wholeupdate = .TRUE. |
---|
546 | if ((deb(n) == -99).AND.(deb(n)==fin(n))) wholeupdate=.TRUE. |
---|
547 | endif |
---|
548 | enddo |
---|
549 | |
---|
550 | IF (present(procname)) THEN |
---|
551 | |
---|
552 | IF (wholeupdate) THEN |
---|
553 | |
---|
554 | Call AGRIF_UpdateWhole |
---|
555 | & (TypeUpdate,parent,child,deb,fin, |
---|
556 | & pttab_Child(1:nbdim),pttab_Parent(1:nbdim), |
---|
557 | & nbtab_Child(1:nbdim),posvartab_Child(1:nbdim), |
---|
558 | & loctab_Child(1:nbdim), |
---|
559 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
560 | & ds_Child(1:nbdim),ds_Parent(1:nbdim),nbdim,procname) |
---|
561 | ELSE |
---|
562 | Call AGRIF_UpdateBcnD |
---|
563 | & (TypeUpdate,parent,child,deb,fin, |
---|
564 | & pttab_Child(1:nbdim),pttab_Parent(1:nbdim), |
---|
565 | & nbtab_Child(1:nbdim),posvartab_Child(1:nbdim), |
---|
566 | & loctab_Child(1:nbdim), |
---|
567 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
568 | & ds_Child(1:nbdim),ds_Parent(1:nbdim),nbdim,procname) |
---|
569 | ENDIF |
---|
570 | ELSE |
---|
571 | IF (wholeupdate) THEN |
---|
572 | Call AGRIF_UpdateWhole |
---|
573 | & (TypeUpdate,parent,child,deb,fin, |
---|
574 | & pttab_Child(1:nbdim),pttab_Parent(1:nbdim), |
---|
575 | & nbtab_Child(1:nbdim),posvartab_Child(1:nbdim), |
---|
576 | & loctab_Child(1:nbdim), |
---|
577 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
578 | & ds_Child(1:nbdim),ds_Parent(1:nbdim),nbdim) |
---|
579 | ELSE |
---|
580 | Call AGRIF_UpdateBcnD |
---|
581 | & (TypeUpdate,parent,child,deb,fin, |
---|
582 | & pttab_Child(1:nbdim),pttab_Parent(1:nbdim), |
---|
583 | & nbtab_Child(1:nbdim),posvartab_Child(1:nbdim), |
---|
584 | & loctab_Child(1:nbdim), |
---|
585 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
586 | & ds_Child(1:nbdim),ds_Parent(1:nbdim),nbdim) |
---|
587 | ENDIF |
---|
588 | ENDIF |
---|
589 | C |
---|
590 | Return |
---|
591 | C |
---|
592 | C |
---|
593 | End subroutine Agrif_UpdateVariable |
---|
594 | C |
---|
595 | C ************************************************************************** |
---|
596 | CCC Subroutine Agrif_UpdateWhole |
---|
597 | C ************************************************************************** |
---|
598 | C |
---|
599 | Subroutine AGRIF_UpdateWhole(TypeUpdate,parent,child,deb,fin, |
---|
600 | & pttab_child,pttab_Parent, |
---|
601 | & nbtab_Child,posvartab_Child, |
---|
602 | & loctab_Child, |
---|
603 | & s_Child,s_Parent, |
---|
604 | & ds_Child,ds_Parent,nbdim,procname) |
---|
605 | C |
---|
606 | CCC Description: |
---|
607 | CCC Subroutine to calculate the boundary conditions for a nD grid variable on |
---|
608 | CCC a fine grid by using a space and time interpolations; it is called by the |
---|
609 | CCC Agrif_CorrectVariable procedure. |
---|
610 | C |
---|
611 | C |
---|
612 | C Declarations: |
---|
613 | C |
---|
614 | |
---|
615 | C |
---|
616 | #ifdef AGRIF_MPI |
---|
617 | C |
---|
618 | #include "mpif.h" |
---|
619 | C |
---|
620 | #endif |
---|
621 | C |
---|
622 | C Arguments |
---|
623 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update (copy or |
---|
624 | ! average) |
---|
625 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent |
---|
626 | ! grid |
---|
627 | TYPE(AGRIF_PVariable) :: child ! Variable on the child |
---|
628 | ! grid |
---|
629 | INTEGER, DIMENSION(6) :: deb, fin |
---|
630 | INTEGER :: nbdim ! Number of dimensions of |
---|
631 | ! the grid variable |
---|
632 | INTEGER,DIMENSION(nbdim) :: pttab_child ! Index of the first point |
---|
633 | ! inside the domain for |
---|
634 | ! the parent grid |
---|
635 | ! variable |
---|
636 | INTEGER,DIMENSION(nbdim) :: pttab_Parent ! Index of the first point |
---|
637 | ! inside the domain for |
---|
638 | ! the child grid |
---|
639 | ! variable |
---|
640 | INTEGER,DIMENSION(nbdim) :: nbtab_Child ! Number of cells of the |
---|
641 | ! child grid |
---|
642 | INTEGER,DIMENSION(nbdim) :: posvartab_Child ! Position of the grid |
---|
643 | ! variable (1 or 2) |
---|
644 | INTEGER,DIMENSION(nbdim) :: loctab_Child ! Indicates if the child |
---|
645 | ! grid has a common |
---|
646 | ! border with the root |
---|
647 | ! grid |
---|
648 | REAL ,DIMENSION(nbdim) :: s_Child,s_Parent ! Positions of the parent |
---|
649 | ! and child grids |
---|
650 | REAL ,DIMENSION(nbdim) :: ds_Child,ds_Parent ! Space steps of the parent |
---|
651 | ! and child grids |
---|
652 | External :: procname |
---|
653 | Optional :: procname |
---|
654 | C |
---|
655 | C Local variables |
---|
656 | INTEGER,DIMENSION(nbdim,2) :: lubglob |
---|
657 | INTEGER :: i |
---|
658 | INTEGER,DIMENSION(nbdim,2,2) :: indtab ! Arrays indicating the |
---|
659 | ! limits of the child |
---|
660 | INTEGER,DIMENSION(nbdim,2,2) :: indtruetab ! grid variable where |
---|
661 | ! boundary conditions are |
---|
662 | integer :: coeffraf |
---|
663 | INTEGER :: debloc, finloc |
---|
664 | C |
---|
665 | #ifdef AGRIF_MPI |
---|
666 | C |
---|
667 | INTEGER,DIMENSION(nbdim) :: lb,ub |
---|
668 | INTEGER,DIMENSION(nbdim,2) :: iminmaxg |
---|
669 | INTEGER :: code |
---|
670 | C |
---|
671 | #endif |
---|
672 | C |
---|
673 | C |
---|
674 | C indtab contains the limits for the fine grid points that will be used |
---|
675 | C in the update scheme |
---|
676 | |
---|
677 | DO i = 1, nbdim |
---|
678 | coeffraf = nint(ds_Parent(i)/ds_Child(i)) |
---|
679 | debloc = 0 |
---|
680 | finloc = nbtab_Child(i)/coeffraf - 1 |
---|
681 | |
---|
682 | IF (posvartab_child(i) == 1) THEN |
---|
683 | finloc = finloc - 1 |
---|
684 | ENDIF |
---|
685 | |
---|
686 | IF (deb(i) > fin(i)) THEN |
---|
687 | debloc = deb(i) |
---|
688 | finloc = finloc - deb(i) |
---|
689 | ENDIF |
---|
690 | |
---|
691 | indtab(i,1,1) = pttab_child(i) + (debloc + 1) * coeffraf |
---|
692 | indtab(i,1,2) = pttab_child(i) + (finloc + 1) * coeffraf |
---|
693 | |
---|
694 | IF (posvartab_child(i) == 1) THEN |
---|
695 | IF (TypeUpdate(i) .EQ. Agrif_Update_Full_Weighting) THEN |
---|
696 | indtab(i,1,1) = indtab(i,1,1) - (coeffraf - 1) |
---|
697 | indtab(i,1,2) = indtab(i,1,2) + (coeffraf - 1) |
---|
698 | ELSE IF (TypeUpdate(i) .NE. Agrif_Update_Copy) THEN |
---|
699 | indtab(i,1,1) = indtab(i,1,1) - coeffraf / 2 |
---|
700 | indtab(i,1,2) = indtab(i,1,2) + coeffraf / 2 |
---|
701 | ENDIF |
---|
702 | ELSE |
---|
703 | indtab(i,1,1) = indtab(i,1,1) - coeffraf |
---|
704 | indtab(i,1,2) = indtab(i,1,2) - 1 |
---|
705 | C at this point, indices are OK for an average |
---|
706 | IF ((TypeUpdate(i) .EQ. Agrif_Update_Full_Weighting)) THEN |
---|
707 | indtab(i,1,1) = indtab(i,1,1) - coeffraf/2 |
---|
708 | indtab(i,1,2) = indtab(i,1,2) + coeffraf/2 |
---|
709 | ENDIF |
---|
710 | ENDIF |
---|
711 | IF (loctab_child(i) == -3) THEN |
---|
712 | indtab(i,1,1) = pttab_child(i) |
---|
713 | C |
---|
714 | if (posvartab_child(i) == 1) then |
---|
715 | C |
---|
716 | indtab(i,1,2) = pttab_child(i) + nbtab_child(i) |
---|
717 | C |
---|
718 | else |
---|
719 | C |
---|
720 | indtab(i,1,2) = pttab_child(i) + nbtab_child(i) - 1 |
---|
721 | ENDIF |
---|
722 | ENDIF |
---|
723 | ENDDO |
---|
724 | |
---|
725 | C lubglob contains the global lbound and ubound of the child array |
---|
726 | C lubglob(:,1) : global lbound for each dimension |
---|
727 | C lubglob(:,2) : global lbound for each dimension |
---|
728 | |
---|
729 | #if !defined AGRIF_MPI |
---|
730 | Call Agrif_nbdim_Get_bound_dimension(child % var,lubglob(:,1), |
---|
731 | & lubglob(:,2),nbdim) |
---|
732 | C |
---|
733 | #else |
---|
734 | C |
---|
735 | Call Agrif_nbdim_Get_bound_dimension(child % var,lb,ub,nbdim) |
---|
736 | DO i = 1,nbdim |
---|
737 | C |
---|
738 | Call Agrif_Invloc(lb(i),Agrif_Procrank,i,iminmaxg(i,1)) |
---|
739 | Call Agrif_Invloc(ub(i),Agrif_Procrank,i,iminmaxg(i,2)) |
---|
740 | C |
---|
741 | ENDDO |
---|
742 | C |
---|
743 | iminmaxg(1:nbdim,2) = - iminmaxg(1:nbdim,2) |
---|
744 | |
---|
745 | CALL MPI_ALLREDUCE(iminmaxg,lubglob,2*nbdim,MPI_INTEGER,MPI_MIN, |
---|
746 | & MPI_COMM_WORLD,code) |
---|
747 | |
---|
748 | lubglob(1:nbdim,2) = - lubglob(1:nbdim,2) |
---|
749 | C |
---|
750 | #endif |
---|
751 | C |
---|
752 | indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1), |
---|
753 | & lubglob(1:nbdim,1)) |
---|
754 | indtruetab(1:nbdim,1,2) = min(indtab(1:nbdim,1,2), |
---|
755 | & lubglob(1:nbdim,2)) |
---|
756 | |
---|
757 | C |
---|
758 | C |
---|
759 | IF (present(procname)) THEN |
---|
760 | Call Agrif_UpdatenD |
---|
761 | & (TypeUpdate,parent,child, |
---|
762 | & indtruetab(1:nbdim,1,1),indtruetab(1:nbdim,1,2), |
---|
763 | & pttab_child(1:nbdim),pttab_Parent(1:nbdim), |
---|
764 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
765 | & ds_Child(1:nbdim),ds_Parent(1:nbdim), |
---|
766 | & posvartab_child,loctab_Child, |
---|
767 | & nbdim,procname) |
---|
768 | ELSE |
---|
769 | Call Agrif_UpdatenD |
---|
770 | & (TypeUpdate,parent,child, |
---|
771 | & indtruetab(1:nbdim,1,1),indtruetab(1:nbdim,1,2), |
---|
772 | & pttab_child(1:nbdim),pttab_Parent(1:nbdim), |
---|
773 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
774 | & ds_Child(1:nbdim),ds_Parent(1:nbdim), |
---|
775 | & posvartab_child,loctab_Child, |
---|
776 | & nbdim) |
---|
777 | ENDIF |
---|
778 | C |
---|
779 | C |
---|
780 | C |
---|
781 | End Subroutine Agrif_UpdateWhole |
---|
782 | C |
---|
783 | C ************************************************************************** |
---|
784 | CCC Subroutine Agrif_UpdateBcnd |
---|
785 | C ************************************************************************** |
---|
786 | C |
---|
787 | Subroutine AGRIF_UpdateBcnd(TypeUpdate,parent,child,deb,fin, |
---|
788 | & pttab_child,pttab_Parent, |
---|
789 | & nbtab_Child,posvartab_Child, |
---|
790 | & loctab_Child, |
---|
791 | & s_Child,s_Parent, |
---|
792 | & ds_Child,ds_Parent,nbdim,procname) |
---|
793 | C |
---|
794 | CCC Description: |
---|
795 | CCC Subroutine to calculate the boundary conditions for a nD grid variable on |
---|
796 | CCC a fine grid by using a space and time interpolations; it is called by the |
---|
797 | CCC Agrif_CorrectVariable procedure. |
---|
798 | C |
---|
799 | C |
---|
800 | C Declarations: |
---|
801 | C |
---|
802 | |
---|
803 | C |
---|
804 | #ifdef AGRIF_MPI |
---|
805 | C |
---|
806 | #include "mpif.h" |
---|
807 | C |
---|
808 | #endif |
---|
809 | C |
---|
810 | C Arguments |
---|
811 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update |
---|
812 | ! (copy or average) |
---|
813 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent |
---|
814 | ! grid |
---|
815 | TYPE(AGRIF_PVariable) :: child ! Variable on the child |
---|
816 | ! grid |
---|
817 | INTEGER, DIMENSION(6) :: deb, fin |
---|
818 | INTEGER :: nbdim ! Number of dimensions of |
---|
819 | ! the grid variable |
---|
820 | INTEGER,DIMENSION(nbdim) :: pttab_child ! Index of the first point |
---|
821 | ! inside the domain for |
---|
822 | ! the parent grid |
---|
823 | ! variable |
---|
824 | INTEGER,DIMENSION(nbdim) :: pttab_Parent ! Index of the first point |
---|
825 | ! inside the domain for |
---|
826 | ! the child grid variable |
---|
827 | INTEGER,DIMENSION(nbdim) :: nbtab_Child ! Number of cells of the |
---|
828 | ! child grid |
---|
829 | INTEGER,DIMENSION(nbdim) :: posvartab_Child ! Position of the grid |
---|
830 | ! variable (1 or 2) |
---|
831 | INTEGER,DIMENSION(nbdim) :: loctab_Child ! Indicates if the child |
---|
832 | ! grid has a common |
---|
833 | ! border with the root |
---|
834 | ! grid |
---|
835 | REAL ,DIMENSION(nbdim) :: s_Child,s_Parent ! Positions of the parent |
---|
836 | ! and child grids |
---|
837 | REAL ,DIMENSION(nbdim) :: ds_Child,ds_Parent ! Space steps of the parent |
---|
838 | ! and child grids |
---|
839 | External :: procname |
---|
840 | Optional :: procname |
---|
841 | C |
---|
842 | C Local variables |
---|
843 | INTEGER,DIMENSION(nbdim,2) :: lubglob |
---|
844 | INTEGER :: i |
---|
845 | INTEGER,DIMENSION(nbdim,2,2) :: indtab ! Arrays indicating the |
---|
846 | ! limits of the child |
---|
847 | INTEGER,DIMENSION(nbdim,2,2) :: indtruetab ! grid variable where |
---|
848 | ! boundary conditions are |
---|
849 | INTEGER,DIMENSION(nbdim,2,2,nbdim) :: ptres ! calculated |
---|
850 | INTEGER :: nb,ndir,n |
---|
851 | integer :: coeffraf |
---|
852 | C |
---|
853 | #ifdef AGRIF_MPI |
---|
854 | C |
---|
855 | INTEGER,DIMENSION(nbdim) :: lb,ub |
---|
856 | INTEGER,DIMENSION(nbdim,2) :: iminmaxg |
---|
857 | INTEGER :: code |
---|
858 | C |
---|
859 | #endif |
---|
860 | C |
---|
861 | C |
---|
862 | |
---|
863 | DO i = 1, nbdim |
---|
864 | coeffraf = nint(ds_Parent(i)/ds_Child(i)) |
---|
865 | indtab(i,1,1) = pttab_child(i) + (deb(i) + 1) * coeffraf |
---|
866 | indtab(i,1,2) = pttab_child(i) + (fin(i) + 1) * coeffraf |
---|
867 | |
---|
868 | indtab(i,2,1) = pttab_child(i) + nbtab_child(i) |
---|
869 | & - (fin(i) + 1) * coeffraf |
---|
870 | indtab(i,2,2) = pttab_child(i) + nbtab_child(i) |
---|
871 | & - (deb(i) + 1) * coeffraf |
---|
872 | |
---|
873 | IF (posvartab_child(i) == 1) THEN |
---|
874 | IF (TypeUpdate(i) .EQ. Agrif_Update_Full_Weighting) THEN |
---|
875 | indtab(i,:,1) = indtab(i,:,1) - (coeffraf - 1) |
---|
876 | indtab(i,:,2) = indtab(i,:,2) + (coeffraf - 1) |
---|
877 | ELSE IF (TypeUpdate(i) .NE. Agrif_Update_Copy) THEN |
---|
878 | indtab(i,:,1) = indtab(i,:,1) - coeffraf / 2 |
---|
879 | indtab(i,:,2) = indtab(i,:,2) + coeffraf / 2 |
---|
880 | ENDIF |
---|
881 | ELSE |
---|
882 | indtab(i,1,1) = indtab(i,1,1) - coeffraf |
---|
883 | indtab(i,1,2) = indtab(i,1,2) - 1 |
---|
884 | indtab(i,2,2) = indtab(i,2,2) + coeffraf - 1 |
---|
885 | IF (TypeUpdate(i) .EQ. Agrif_Update_Full_Weighting) THEN |
---|
886 | indtab(i,1,1) = indtab(i,1,1) - coeffraf/2 |
---|
887 | indtab(i,1,2) = indtab(i,1,2) + coeffraf/2 |
---|
888 | indtab(i,2,1) = indtab(i,2,1) - coeffraf/2 |
---|
889 | indtab(i,2,2) = indtab(i,2,2) + coeffraf/2 |
---|
890 | ENDIF |
---|
891 | ENDIF |
---|
892 | ENDDO |
---|
893 | |
---|
894 | #if !defined AGRIF_MPI |
---|
895 | Call Agrif_nbdim_Get_bound_dimension(child % var,lubglob(:,1), |
---|
896 | & lubglob(:,2),nbdim) |
---|
897 | |
---|
898 | C |
---|
899 | #else |
---|
900 | C |
---|
901 | Call Agrif_nbdim_Get_bound_dimension(child % var,lb,ub,nbdim) |
---|
902 | DO i = 1,nbdim |
---|
903 | C |
---|
904 | Call Agrif_Invloc(lb(i),Agrif_Procrank,i,iminmaxg(i,1)) |
---|
905 | Call Agrif_Invloc(ub(i),Agrif_Procrank,i,iminmaxg(i,2)) |
---|
906 | C |
---|
907 | ENDDO |
---|
908 | C |
---|
909 | iminmaxg(1:nbdim,2) = - iminmaxg(1:nbdim,2) |
---|
910 | |
---|
911 | CALL MPI_ALLREDUCE(iminmaxg,lubglob,2*nbdim,MPI_INTEGER,MPI_MIN, |
---|
912 | & MPI_COMM_WORLD,code) |
---|
913 | |
---|
914 | lubglob(1:nbdim,2) = - lubglob(1:nbdim,2) |
---|
915 | C |
---|
916 | #endif |
---|
917 | C |
---|
918 | indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1), |
---|
919 | & lubglob(1:nbdim,1)) |
---|
920 | indtruetab(1:nbdim,1,2) = max(indtab(1:nbdim,1,2), |
---|
921 | & lubglob(1:nbdim,1)) |
---|
922 | indtruetab(1:nbdim,2,1) = min(indtab(1:nbdim,2,1), |
---|
923 | & lubglob(1:nbdim,2)) |
---|
924 | indtruetab(1:nbdim,2,2) = min(indtab(1:nbdim,2,2), |
---|
925 | & lubglob(1:nbdim,2)) |
---|
926 | |
---|
927 | C |
---|
928 | C |
---|
929 | do nb = 1,nbdim |
---|
930 | C |
---|
931 | do ndir = 1,2 |
---|
932 | C |
---|
933 | if (loctab_child(nb) /= -3) then |
---|
934 | C |
---|
935 | do n = 1,2 |
---|
936 | C |
---|
937 | ptres(nb,n,ndir,nb) = indtruetab(nb,ndir,n) |
---|
938 | C |
---|
939 | enddo |
---|
940 | C |
---|
941 | do i = 1,nbdim |
---|
942 | C |
---|
943 | if (i .NE. nb) then |
---|
944 | C |
---|
945 | if (loctab_child(i) == -3) then |
---|
946 | C |
---|
947 | ptres(i,1,ndir,nb) = pttab_child(i) |
---|
948 | C |
---|
949 | else |
---|
950 | C |
---|
951 | ptres(i,1,ndir,nb) = indtruetab(i,1,1) |
---|
952 | C |
---|
953 | endif |
---|
954 | C |
---|
955 | if (loctab_child(i) == -3) then |
---|
956 | C |
---|
957 | if (posvartab_child(i) == 1) then |
---|
958 | C |
---|
959 | ptres(i,2,ndir,nb) = pttab_child(i) |
---|
960 | & + nbtab_child(i) |
---|
961 | C |
---|
962 | else |
---|
963 | C |
---|
964 | ptres(i,2,ndir,nb) = pttab_child(i) |
---|
965 | & + nbtab_child(i) - 1 |
---|
966 | C |
---|
967 | endif |
---|
968 | C |
---|
969 | else |
---|
970 | C |
---|
971 | ptres(i,2,ndir,nb) = indtruetab(i,2,2) |
---|
972 | C |
---|
973 | endif |
---|
974 | C |
---|
975 | endif |
---|
976 | C |
---|
977 | enddo |
---|
978 | |
---|
979 | C |
---|
980 | |
---|
981 | endif |
---|
982 | |
---|
983 | enddo |
---|
984 | enddo |
---|
985 | C |
---|
986 | |
---|
987 | C |
---|
988 | |
---|
989 | do nb = 1,nbdim |
---|
990 | C |
---|
991 | do ndir = 1,2 |
---|
992 | C |
---|
993 | if (loctab_child(nb) /= -3) then |
---|
994 | C |
---|
995 | IF (present(procname)) THEN |
---|
996 | Call Agrif_UpdatenD |
---|
997 | & (TypeUpdate,parent,child, |
---|
998 | & ptres(1:nbdim,1,ndir,nb),ptres(1:nbdim,2,ndir,nb), |
---|
999 | & pttab_child(1:nbdim),pttab_Parent(1:nbdim), |
---|
1000 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
1001 | & ds_Child(1:nbdim),ds_Parent(1:nbdim), |
---|
1002 | & posvartab_Child,loctab_Child, |
---|
1003 | & nbdim,procname,nb,ndir) |
---|
1004 | ELSE |
---|
1005 | Call Agrif_UpdatenD |
---|
1006 | & (TypeUpdate,parent,child, |
---|
1007 | & ptres(1:nbdim,1,ndir,nb),ptres(1:nbdim,2,ndir,nb), |
---|
1008 | & pttab_child(1:nbdim),pttab_Parent(1:nbdim), |
---|
1009 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
1010 | & ds_Child(1:nbdim),ds_Parent(1:nbdim), |
---|
1011 | & posvartab_Child,loctab_Child, |
---|
1012 | & nbdim) |
---|
1013 | ENDIF |
---|
1014 | C |
---|
1015 | endif |
---|
1016 | |
---|
1017 | C |
---|
1018 | enddo |
---|
1019 | C |
---|
1020 | enddo |
---|
1021 | C |
---|
1022 | C |
---|
1023 | C |
---|
1024 | End Subroutine Agrif_UpdateBcnd |
---|
1025 | C |
---|
1026 | C ************************************************************************** |
---|
1027 | CCC Subroutine Agrif_UpdatenD |
---|
1028 | C ************************************************************************** |
---|
1029 | C |
---|
1030 | Subroutine Agrif_UpdatenD(TypeUpdate,parent,child, |
---|
1031 | & pttab,petab, |
---|
1032 | & pttab_Child,pttab_Parent, |
---|
1033 | & s_Child,s_Parent, |
---|
1034 | & ds_Child,ds_Parent, |
---|
1035 | & posvartab_Child,loctab_Child, |
---|
1036 | & nbdim,procname,nb,ndir) |
---|
1037 | C |
---|
1038 | C Description: |
---|
1039 | C Subroutine to update a 2D grid variable on the parent grid of |
---|
1040 | C the current grid. |
---|
1041 | C |
---|
1042 | C Declarations: |
---|
1043 | C |
---|
1044 | |
---|
1045 | C |
---|
1046 | #ifdef AGRIF_MPI |
---|
1047 | C |
---|
1048 | #include "mpif.h" |
---|
1049 | C |
---|
1050 | #endif |
---|
1051 | C |
---|
1052 | C Arguments |
---|
1053 | INTEGER :: nbdim |
---|
1054 | INTEGER, DIMENSION(6) :: TypeUpdate ! TYPE of update |
---|
1055 | ! (copy or average) |
---|
1056 | TYPE(AGRIF_PVARIABLE) :: parent ! Variable of the parent |
---|
1057 | ! grid |
---|
1058 | TYPE(AGRIF_PVARIABLE) :: child ! Variable of the child |
---|
1059 | ! grid |
---|
1060 | INTEGER,DIMENSION(nbdim) :: pttab ! Index of the first |
---|
1061 | ! point inside the |
---|
1062 | ! domain |
---|
1063 | INTEGER,DIMENSION(nbdim) :: petab ! Index of the first |
---|
1064 | ! point inside the |
---|
1065 | ! domain |
---|
1066 | INTEGER,DIMENSION(nbdim) :: pttab_Child ! Index of the first |
---|
1067 | ! point inside the |
---|
1068 | ! domain for the child |
---|
1069 | ! grid variable |
---|
1070 | INTEGER,DIMENSION(nbdim) :: pttab_Parent ! Index of the first |
---|
1071 | ! point inside the |
---|
1072 | ! domain for the parent |
---|
1073 | ! grid variable |
---|
1074 | REAL,DIMENSION(nbdim) :: s_Child,s_Parent ! Positions of the parent |
---|
1075 | ! and child grids |
---|
1076 | REAL,DIMENSION(nbdim) :: ds_Child,ds_Parent ! Space steps of the |
---|
1077 | ! parent and child |
---|
1078 | ! grids |
---|
1079 | External :: procname |
---|
1080 | Optional :: procname |
---|
1081 | Integer :: nb,ndir |
---|
1082 | Optional :: nb,ndir |
---|
1083 | |
---|
1084 | C |
---|
1085 | C Local pointers |
---|
1086 | TYPE(AGRIF_PVARIABLE), SAVE :: tempP ! Temporary parent grid variable |
---|
1087 | TYPE(AGRIF_PVARIABLE), SAVE :: tempC ! Temporary child grid variable |
---|
1088 | C |
---|
1089 | C Local scalars |
---|
1090 | INTEGER,DIMENSION(nbdim) :: pttruetab,cetruetab |
---|
1091 | INTEGER,DIMENSION(nbdim) :: posvartab_Child,loctab_Child |
---|
1092 | INTEGER,DIMENSION(nbdim) :: indmin,indmax |
---|
1093 | INTEGER,DIMENSION(nbdim) :: indminglob,indmaxglob |
---|
1094 | REAL ,DIMENSION(nbdim) :: s_Child_temp,s_Parent_temp |
---|
1095 | cccccccc LOGICAL,DIMENSION(nbdim) :: noraftab |
---|
1096 | INTEGER,DIMENSION(nbdim) :: lowerbound,upperbound |
---|
1097 | LOGICAL :: memberin, member |
---|
1098 | INTEGER,DIMENSION(nbdim) :: pttruetabwhole,cetruetabwhole |
---|
1099 | INTEGER,DIMENSION(nbdim,2,2) :: childarray |
---|
1100 | INTEGER,DIMENSION(nbdim,2,2) :: parentarray |
---|
1101 | TYPE(AGRIF_PVARIABLE), SAVE :: tempCextend,tempPextend ! Temporary child |
---|
1102 | INTEGER :: nbin, ndirin |
---|
1103 | C |
---|
1104 | #ifdef AGRIF_MPI |
---|
1105 | C |
---|
1106 | INTEGER,DIMENSION(nbdim) :: indminglob2,indmaxglob2 |
---|
1107 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: sendtoproc1,recvfromproc1 |
---|
1108 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: sendtoproc2,recvfromproc2 |
---|
1109 | INTEGER :: code |
---|
1110 | INTEGER :: i,j,k |
---|
1111 | INTEGER,DIMENSION(nbdim,4) :: tab3 |
---|
1112 | INTEGER,DIMENSION(nbdim,4,0:Agrif_Nbprocs-1) :: tab4 |
---|
1113 | INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,8) :: tab4t |
---|
1114 | INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,8) :: tab5t |
---|
1115 | LOGICAL :: find_list_update |
---|
1116 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: memberinall, memberinall2 |
---|
1117 | LOGICAL, DIMENSION(1) :: memberin1 |
---|
1118 | C |
---|
1119 | #endif |
---|
1120 | C |
---|
1121 | |
---|
1122 | C |
---|
1123 | C local lbound and ubound of the child array |
---|
1124 | |
---|
1125 | Call Agrif_nbdim_Get_bound_dimension(child%var, |
---|
1126 | & lowerbound,upperbound,nbdim) |
---|
1127 | |
---|
1128 | C here pttab and petab corresponds to the (global) indices of the points needed |
---|
1129 | C in the update |
---|
1130 | C pttruetab and cetruetab contains only indices that are present |
---|
1131 | C on the local processor |
---|
1132 | |
---|
1133 | Call Agrif_Childbounds(nbdim, |
---|
1134 | & lowerbound,upperbound, |
---|
1135 | & pttab,petab, |
---|
1136 | & pttruetab,cetruetab,memberin) |
---|
1137 | |
---|
1138 | Call Agrif_Prtbounds(nbdim,indminglob,indmaxglob,s_Parent_temp, |
---|
1139 | & s_Child_temp,s_Child,ds_Child, |
---|
1140 | & s_Parent,ds_Parent, |
---|
1141 | & pttab,petab,pttab_Child, |
---|
1142 | & pttab_Parent, |
---|
1143 | & posvartab_Child,TypeUpdate,loctab_Child |
---|
1144 | #ifdef AGRIF_MPI |
---|
1145 | & ,pttruetabwhole,cetruetabwhole |
---|
1146 | #endif |
---|
1147 | & ) |
---|
1148 | |
---|
1149 | #ifdef AGRIF_MPI |
---|
1150 | IF (memberin) THEN |
---|
1151 | Call Agrif_GlobtoLocInd2(childarray, |
---|
1152 | & lowerbound,upperbound, |
---|
1153 | & pttruetab,cetruetab, |
---|
1154 | & nbdim,Agrif_Procrank, |
---|
1155 | & member) |
---|
1156 | |
---|
1157 | ENDIF |
---|
1158 | |
---|
1159 | |
---|
1160 | Call Agrif_Prtbounds(nbdim,indmin,indmax,s_Parent_temp, |
---|
1161 | & s_Child_temp,s_Child,ds_Child, |
---|
1162 | & s_Parent,ds_Parent, |
---|
1163 | & pttruetab,cetruetab,pttab_Child, |
---|
1164 | & pttab_Parent, |
---|
1165 | & posvartab_Child,TypeUpdate,loctab_Child |
---|
1166 | & ,pttruetabwhole,cetruetabwhole |
---|
1167 | & ) |
---|
1168 | |
---|
1169 | #else |
---|
1170 | indmin = indminglob |
---|
1171 | indmax = indmaxglob |
---|
1172 | pttruetabwhole = pttruetab |
---|
1173 | cetruetabwhole = cetruetab |
---|
1174 | childarray(:,1,2) = pttruetab |
---|
1175 | childarray(:,2,2) = cetruetab |
---|
1176 | #endif |
---|
1177 | |
---|
1178 | IF (present(procname)) THEN |
---|
1179 | IF (.Not.present(nb)) THEN |
---|
1180 | nbin=0 |
---|
1181 | ndirin=0 |
---|
1182 | ELSE |
---|
1183 | nbin = nb |
---|
1184 | ndirin = ndir |
---|
1185 | ENDIF |
---|
1186 | ENDIF |
---|
1187 | |
---|
1188 | IF (memberin) THEN |
---|
1189 | IF (.not.associated(tempC%var)) allocate(tempC%var) |
---|
1190 | |
---|
1191 | C |
---|
1192 | Call Agrif_nbdim_allocation(tempC%var, |
---|
1193 | & pttruetab,cetruetab,nbdim) |
---|
1194 | |
---|
1195 | Call Agrif_nbdim_Full_VarEQreal(tempC%var,0.,nbdim) |
---|
1196 | |
---|
1197 | IF (present(procname)) THEN |
---|
1198 | SELECT CASE (nbdim) |
---|
1199 | CASE(1) |
---|
1200 | CALL procname(tempC%var%array1, |
---|
1201 | & childarray(1,1,2),childarray(1,2,2), |
---|
1202 | & .TRUE.,nbin,ndirin) |
---|
1203 | CASE(2) |
---|
1204 | CALL procname(tempC%var%array2, |
---|
1205 | & childarray(1,1,2),childarray(1,2,2), |
---|
1206 | & childarray(2,1,2),childarray(2,2,2), |
---|
1207 | & .TRUE.,nbin,ndirin) |
---|
1208 | CASE(3) |
---|
1209 | CALL procname(tempC%var%array3, |
---|
1210 | & childarray(1,1,2),childarray(1,2,2), |
---|
1211 | & childarray(2,1,2),childarray(2,2,2), |
---|
1212 | & childarray(3,1,2),childarray(3,2,2), |
---|
1213 | & .TRUE.,nbin,ndirin) |
---|
1214 | CASE(4) |
---|
1215 | CALL procname(tempC%var%array4, |
---|
1216 | & 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 | & childarray(4,1,2),childarray(4,2,2), |
---|
1220 | & .TRUE.,nbin,ndirin) |
---|
1221 | CASE(5) |
---|
1222 | CALL procname(tempC%var%array5, |
---|
1223 | & childarray(1,1,2),childarray(1,2,2), |
---|
1224 | & childarray(2,1,2),childarray(2,2,2), |
---|
1225 | & childarray(3,1,2),childarray(3,2,2), |
---|
1226 | & childarray(4,1,2),childarray(4,2,2), |
---|
1227 | & childarray(5,1,2),childarray(5,2,2), |
---|
1228 | & .TRUE.,nbin,ndirin) |
---|
1229 | CASE(6) |
---|
1230 | CALL procname(tempC%var%array6, |
---|
1231 | & childarray(1,1,2),childarray(1,2,2), |
---|
1232 | & childarray(2,1,2),childarray(2,2,2), |
---|
1233 | & childarray(3,1,2),childarray(3,2,2), |
---|
1234 | & childarray(4,1,2),childarray(4,2,2), |
---|
1235 | & childarray(5,1,2),childarray(5,2,2), |
---|
1236 | & childarray(6,1,2),childarray(6,2,2), |
---|
1237 | & .TRUE.,nbin,ndirin) |
---|
1238 | END SELECT |
---|
1239 | ELSE |
---|
1240 | Call Agrif_nbdim_VarEQvar(tempC%var,pttruetab,cetruetab, |
---|
1241 | & child%var,childarray(:,1,2),childarray(:,2,2), |
---|
1242 | & nbdim) |
---|
1243 | ENDIF |
---|
1244 | |
---|
1245 | ENDIF |
---|
1246 | |
---|
1247 | |
---|
1248 | |
---|
1249 | C |
---|
1250 | C |
---|
1251 | #ifdef AGRIF_MPI |
---|
1252 | C |
---|
1253 | C tab2 contains the necessary limits of the parent grid for each processor |
---|
1254 | |
---|
1255 | IF (Associated(child%var%list_update)) THEN |
---|
1256 | Call Agrif_Find_list_update(child%var%list_update,pttab,petab, |
---|
1257 | & pttab_Child,pttab_Parent,nbdim, |
---|
1258 | & find_list_update,tab4t,tab5t,memberinall,memberinall2, |
---|
1259 | & sendtoproc1,recvfromproc1,sendtoproc2,recvfromproc2) |
---|
1260 | ELSE |
---|
1261 | find_list_update = .FALSE. |
---|
1262 | ENDIF |
---|
1263 | |
---|
1264 | if (.not.find_list_update) then |
---|
1265 | tab3(:,1) = pttruetab(:) |
---|
1266 | tab3(:,2) = cetruetab(:) |
---|
1267 | tab3(:,3) = pttruetabwhole(:) |
---|
1268 | tab3(:,4) = cetruetabwhole(:) |
---|
1269 | C |
---|
1270 | C |
---|
1271 | Call MPI_ALLGATHER(tab3,4*nbdim,MPI_INTEGER,tab4,4*nbdim, |
---|
1272 | & MPI_INTEGER,MPI_COMM_WORLD,code) |
---|
1273 | |
---|
1274 | IF (.not.associated(tempCextend%var)) Allocate(tempCextend%var) |
---|
1275 | DO k=0,Agrif_Nbprocs-1 |
---|
1276 | do j=1,4 |
---|
1277 | do i=1,nbdim |
---|
1278 | tab4t(i,k,j) = tab4(i,j,k) |
---|
1279 | enddo |
---|
1280 | enddo |
---|
1281 | enddo |
---|
1282 | |
---|
1283 | memberin1(1) = memberin |
---|
1284 | CALL MPI_ALLGATHER(memberin1,1,MPI_LOGICAL,memberinall, |
---|
1285 | & 1,MPI_LOGICAL,MPI_COMM_WORLD,code) |
---|
1286 | |
---|
1287 | Call Get_External_Data_first(tab4t(:,:,1), |
---|
1288 | & tab4t(:,:,2), |
---|
1289 | & tab4t(:,:,3),tab4t(:,:,4),nbdim,memberin,memberin, |
---|
1290 | & memberinall,sendtoproc1,recvfromproc1,tab4t(:,:,5), |
---|
1291 | & tab4t(:,:,6),tab4t(:,:,7),tab4t(:,:,8)) |
---|
1292 | |
---|
1293 | endif |
---|
1294 | |
---|
1295 | Call ExchangeSameLevel2(sendtoproc1,recvfromproc1,nbdim, |
---|
1296 | & tab4t(:,:,3),tab4t(:,:,4),tab4t(:,:,5),tab4t(:,:,6), |
---|
1297 | & tab4t(:,:,7),tab4t(:,:,8),memberin,tempC, |
---|
1298 | & tempCextend) |
---|
1299 | |
---|
1300 | ! Call Get_External_Data(tempC,tempCextend,tab4t(:,:,1), |
---|
1301 | ! & tab4t(:,:,2), |
---|
1302 | ! & tab4t(:,:,3),tab4t(:,:,4),nbdim,memberin,memberin, |
---|
1303 | ! & memberinall) |
---|
1304 | |
---|
1305 | #else |
---|
1306 | tempCextend%var => tempC%var |
---|
1307 | #endif |
---|
1308 | |
---|
1309 | C |
---|
1310 | C |
---|
1311 | C Update of the parent grid (tempP) from the child grid (tempC) |
---|
1312 | |
---|
1313 | |
---|
1314 | IF (memberin) THEN |
---|
1315 | |
---|
1316 | IF (.not.associated(tempP%var)) allocate(tempP%var) |
---|
1317 | Call Agrif_nbdim_allocation(tempP%var, |
---|
1318 | & indmin,indmax,nbdim) |
---|
1319 | |
---|
1320 | if ( nbdim .EQ. 1 ) then |
---|
1321 | tempP%var%array1 = 0. |
---|
1322 | Call Agrif_Update_1D_recursive(TypeUpdate, |
---|
1323 | & tempP%var%array1,tempCextend%var%array1, |
---|
1324 | & indmin,indmax, |
---|
1325 | & pttruetabwhole,cetruetabwhole, |
---|
1326 | & s_Child_temp,s_Parent_temp, |
---|
1327 | & ds_Child,ds_Parent,nbdim) |
---|
1328 | endif |
---|
1329 | if ( nbdim .EQ. 2 ) then |
---|
1330 | Call Agrif_Update_2D_recursive(TypeUpdate, |
---|
1331 | & tempP%var%array2,tempCextend%var%array2, |
---|
1332 | & indmin,indmax, |
---|
1333 | & pttruetabwhole,cetruetabwhole, |
---|
1334 | & s_Child_temp,s_Parent_temp, |
---|
1335 | & ds_Child,ds_Parent,nbdim) |
---|
1336 | endif |
---|
1337 | |
---|
1338 | if ( nbdim .EQ. 3 ) then |
---|
1339 | Call Agrif_Update_3D_recursive(TypeUpdate, |
---|
1340 | & tempP%var%array3,tempCextend%var%array3, |
---|
1341 | & indmin,indmax, |
---|
1342 | & pttruetabwhole,cetruetabwhole, |
---|
1343 | & s_Child_temp,s_Parent_temp, |
---|
1344 | & ds_Child,ds_Parent,nbdim) |
---|
1345 | endif |
---|
1346 | if ( nbdim .EQ. 4 ) then |
---|
1347 | Call Agrif_Update_4D_recursive(TypeUpdate, |
---|
1348 | & tempP%var%array4,tempCextend%var%array4, |
---|
1349 | & indmin,indmax, |
---|
1350 | & pttruetabwhole,cetruetabwhole, |
---|
1351 | & s_Child_temp,s_Parent_temp, |
---|
1352 | & ds_Child,ds_Parent,nbdim) |
---|
1353 | endif |
---|
1354 | if ( nbdim .EQ. 5 ) then |
---|
1355 | Call Agrif_Update_5D_recursive(TypeUpdate, |
---|
1356 | & tempP%var%array5,tempCextend%var%array5, |
---|
1357 | & indmin,indmax, |
---|
1358 | & pttruetabwhole,cetruetabwhole, |
---|
1359 | & s_Child_temp,s_Parent_temp, |
---|
1360 | & ds_Child,ds_Parent,nbdim) |
---|
1361 | endif |
---|
1362 | if ( nbdim .EQ. 6 ) then |
---|
1363 | Call Agrif_Update_6D_recursive(TypeUpdate, |
---|
1364 | & tempP%var%array6,tempCextend%var%array6, |
---|
1365 | & indmin,indmax, |
---|
1366 | & pttruetabwhole,cetruetabwhole, |
---|
1367 | & s_Child_temp,s_Parent_temp, |
---|
1368 | & ds_Child,ds_Parent,nbdim) |
---|
1369 | endif |
---|
1370 | |
---|
1371 | Call Agrif_nbdim_deallocation(tempCextend%var,nbdim) |
---|
1372 | C Deallocate(tempCextend%var) |
---|
1373 | |
---|
1374 | ENDIF |
---|
1375 | |
---|
1376 | #ifdef AGRIF_MPI |
---|
1377 | Call Agrif_nbdim_Get_bound_dimension(parent%var, |
---|
1378 | & lowerbound,upperbound,nbdim) |
---|
1379 | |
---|
1380 | Call Agrif_ChildGrid_to_ParentGrid() |
---|
1381 | C |
---|
1382 | Call Agrif_Childbounds(nbdim, |
---|
1383 | & lowerbound,upperbound, |
---|
1384 | & indminglob,indmaxglob, |
---|
1385 | & indminglob2,indmaxglob2,member) |
---|
1386 | C |
---|
1387 | IF (member) THEN |
---|
1388 | Call Agrif_GlobtoLocInd2(parentarray, |
---|
1389 | & lowerbound,upperbound, |
---|
1390 | & indminglob2,indmaxglob2, |
---|
1391 | & nbdim,Agrif_Procrank, |
---|
1392 | & member) |
---|
1393 | ENDIF |
---|
1394 | |
---|
1395 | Call Agrif_ParentGrid_to_ChildGrid() |
---|
1396 | |
---|
1397 | if (.not.find_list_update) then |
---|
1398 | tab3(:,1) = indmin(:) |
---|
1399 | tab3(:,2) = indmax(:) |
---|
1400 | tab3(:,3) = indminglob2(:) |
---|
1401 | tab3(:,4) = indmaxglob2(:) |
---|
1402 | C |
---|
1403 | Call MPI_ALLGATHER(tab3,4*nbdim,MPI_INTEGER,tab4,4*nbdim, |
---|
1404 | & MPI_INTEGER,MPI_COMM_WORLD,code) |
---|
1405 | |
---|
1406 | IF (.not.associated(tempPextend%var)) Allocate(tempPextend%var) |
---|
1407 | DO k=0,Agrif_Nbprocs-1 |
---|
1408 | do j=1,4 |
---|
1409 | do i=1,nbdim |
---|
1410 | tab5t(i,k,j) = tab4(i,j,k) |
---|
1411 | enddo |
---|
1412 | enddo |
---|
1413 | enddo |
---|
1414 | |
---|
1415 | memberin1(1) = member |
---|
1416 | CALL MPI_ALLGATHER(memberin1,1,MPI_LOGICAL,memberinall2, |
---|
1417 | & 1,MPI_LOGICAL,MPI_COMM_WORLD,code) |
---|
1418 | |
---|
1419 | Call Get_External_Data_first(tab5t(:,:,1), |
---|
1420 | & tab5t(:,:,2), |
---|
1421 | & tab5t(:,:,3),tab5t(:,:,4),nbdim,memberin,member, |
---|
1422 | & memberinall2,sendtoproc2,recvfromproc2,tab5t(:,:,5), |
---|
1423 | & tab5t(:,:,6),tab5t(:,:,7),tab5t(:,:,8)) |
---|
1424 | |
---|
1425 | Call Agrif_Addto_list_update(child%var%list_update,pttab,petab, |
---|
1426 | & pttab_Child,pttab_Parent,nbdim |
---|
1427 | & ,tab4t,tab5t,memberinall,memberinall2, |
---|
1428 | & sendtoproc1,recvfromproc1,sendtoproc2,recvfromproc2) |
---|
1429 | |
---|
1430 | endif |
---|
1431 | |
---|
1432 | c Call Get_External_Data(tempP,tempPextend,tab5t(:,:,1), |
---|
1433 | c & tab5t(:,:,2), |
---|
1434 | c & tab5t(:,:,3),tab5t(:,:,4),nbdim,memberin,member, |
---|
1435 | c & memberinall2) |
---|
1436 | |
---|
1437 | Call ExchangeSameLevel2(sendtoproc2,recvfromproc2,nbdim, |
---|
1438 | & tab5t(:,:,3),tab5t(:,:,4),tab5t(:,:,5),tab5t(:,:,6), |
---|
1439 | & tab5t(:,:,7),tab5t(:,:,8),member,tempP, |
---|
1440 | & tempPextend) |
---|
1441 | |
---|
1442 | #else |
---|
1443 | tempPextend%var => tempP%var |
---|
1444 | parentarray(:,1,1) = indmin |
---|
1445 | parentarray(:,2,1) = indmax |
---|
1446 | parentarray(:,1,2) = indmin |
---|
1447 | parentarray(:,2,2) = indmax |
---|
1448 | member = .TRUE. |
---|
1449 | #endif |
---|
1450 | |
---|
1451 | C |
---|
1452 | C |
---|
1453 | C |
---|
1454 | C Special values on the child grid |
---|
1455 | if (Agrif_UseSpecialValueFineGrid) then |
---|
1456 | C |
---|
1457 | ccc noraftab(1:nbdim) = |
---|
1458 | ccc & child % var % root_var % interptab(1:nbdim) .EQ. 'N' |
---|
1459 | C |
---|
1460 | #ifdef AGRIF_MPI |
---|
1461 | C |
---|
1462 | c Allocate(childvalues% var) |
---|
1463 | C |
---|
1464 | c Call Agrif_nbdim_allocation(childvalues%var, |
---|
1465 | c & pttruetab,cetruetab,nbdim) |
---|
1466 | c Call Agrif_nbdim_Full_VarEQvar(childvalues% var, |
---|
1467 | c & tempC% var, |
---|
1468 | c & nbdim) |
---|
1469 | c Call Agrif_CheckMasknD(tempC,childvalues, |
---|
1470 | c & pttruetab(1:nbdim),cetruetab(1:nbdim), |
---|
1471 | c & pttruetab(1:nbdim),cetruetab(1:nbdim), |
---|
1472 | c & noraftab(1:nbdim),nbdim) |
---|
1473 | c Call Agrif_nbdim_deallocation(childvalues% var,nbdim) |
---|
1474 | c Deallocate(childvalues % var) |
---|
1475 | C |
---|
1476 | #else |
---|
1477 | C |
---|
1478 | c Call Agrif_nbdim_Get_bound_dimension(child%var, |
---|
1479 | c & lowerbound,upperbound,nbdim) |
---|
1480 | c Call Agrif_CheckMasknD(tempC,child, |
---|
1481 | c & pttruetab(1:nbdim),cetruetab(1:nbdim), |
---|
1482 | c & lowerbound, |
---|
1483 | c & upperbound, |
---|
1484 | c & noraftab(1:nbdim),nbdim) |
---|
1485 | C |
---|
1486 | #endif |
---|
1487 | C |
---|
1488 | endif |
---|
1489 | |
---|
1490 | |
---|
1491 | C |
---|
1492 | C |
---|
1493 | C |
---|
1494 | C |
---|
1495 | C Special values on the parent grid |
---|
1496 | if (Agrif_UseSpecialValue) then |
---|
1497 | C |
---|
1498 | #ifdef AGRIF_MPI |
---|
1499 | C |
---|
1500 | c Call GiveAgrif_SpecialValueToTab_mpi(parent%var,tempP%var, |
---|
1501 | c & parentarray, |
---|
1502 | c & indmin,indmax, |
---|
1503 | c & Agrif_SpecialValue,nbdim) |
---|
1504 | C |
---|
1505 | C |
---|
1506 | #else |
---|
1507 | C |
---|
1508 | c Call GiveAgrif_SpecialValueToTab(parent%var,tempP%var, |
---|
1509 | c & indmin,indmax, |
---|
1510 | c & Agrif_SpecialValue,nbdim) |
---|
1511 | C |
---|
1512 | #endif |
---|
1513 | C |
---|
1514 | C |
---|
1515 | endif |
---|
1516 | C |
---|
1517 | C |
---|
1518 | IF (member) THEN |
---|
1519 | |
---|
1520 | IF (present(procname)) THEN |
---|
1521 | CALL Agrif_ChildGrid_to_ParentGrid() |
---|
1522 | SELECT CASE(nbdim) |
---|
1523 | CASE(1) |
---|
1524 | CALL procname( |
---|
1525 | & tempPextend%var%array1( |
---|
1526 | & parentarray(1,1,1):parentarray(1,2,1)), |
---|
1527 | & parentarray(1,1,2),parentarray(1,2,2), |
---|
1528 | & .FALSE.,nbin,ndirin |
---|
1529 | & ) |
---|
1530 | CASE(2) |
---|
1531 | CALL procname( |
---|
1532 | & tempPextend%var%array2( |
---|
1533 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1534 | & parentarray(2,1,1):parentarray(2,2,1)), |
---|
1535 | & parentarray(1,1,2),parentarray(1,2,2), |
---|
1536 | & parentarray(2,1,2),parentarray(2,2,2), |
---|
1537 | & .FALSE.,nbin,ndirin |
---|
1538 | & ) |
---|
1539 | CASE(3) |
---|
1540 | CALL procname( |
---|
1541 | & tempPextend%var%array3( |
---|
1542 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1543 | & parentarray(2,1,1):parentarray(2,2,1), |
---|
1544 | & parentarray(3,1,1):parentarray(3,2,1)), |
---|
1545 | & parentarray(1,1,2),parentarray(1,2,2), |
---|
1546 | & parentarray(2,1,2),parentarray(2,2,2), |
---|
1547 | & parentarray(3,1,2),parentarray(3,2,2), |
---|
1548 | & .FALSE.,nbin,ndirin |
---|
1549 | & ) |
---|
1550 | CASE(4) |
---|
1551 | CALL procname( |
---|
1552 | & tempPextend%var%array4( |
---|
1553 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1554 | & parentarray(2,1,1):parentarray(2,2,1), |
---|
1555 | & parentarray(3,1,1):parentarray(3,2,1), |
---|
1556 | & parentarray(4,1,1):parentarray(4,2,1)), |
---|
1557 | & parentarray(1,1,2),parentarray(1,2,2), |
---|
1558 | & parentarray(2,1,2),parentarray(2,2,2), |
---|
1559 | & parentarray(3,1,2),parentarray(3,2,2), |
---|
1560 | & parentarray(4,1,2),parentarray(4,2,2), |
---|
1561 | & .FALSE.,nbin,ndirin |
---|
1562 | & ) |
---|
1563 | CASE(5) |
---|
1564 | CALL procname( |
---|
1565 | & tempPextend%var%array5( |
---|
1566 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1567 | & parentarray(2,1,1):parentarray(2,2,1), |
---|
1568 | & parentarray(3,1,1):parentarray(3,2,1), |
---|
1569 | & parentarray(4,1,1):parentarray(4,2,1), |
---|
1570 | & parentarray(5,1,1):parentarray(5,2,1)), |
---|
1571 | & parentarray(1,1,2),parentarray(1,2,2), |
---|
1572 | & parentarray(2,1,2),parentarray(2,2,2), |
---|
1573 | & parentarray(3,1,2),parentarray(3,2,2), |
---|
1574 | & parentarray(4,1,2),parentarray(4,2,2), |
---|
1575 | & parentarray(5,1,2),parentarray(5,2,2), |
---|
1576 | & .FALSE.,nbin,ndirin |
---|
1577 | & ) |
---|
1578 | CASE(6) |
---|
1579 | CALL procname( |
---|
1580 | & tempPextend%var%array6( |
---|
1581 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1582 | & parentarray(2,1,1):parentarray(2,2,1), |
---|
1583 | & parentarray(3,1,1):parentarray(3,2,1), |
---|
1584 | & parentarray(4,1,1):parentarray(4,2,1), |
---|
1585 | & parentarray(5,1,1):parentarray(5,2,1), |
---|
1586 | & parentarray(6,1,1):parentarray(6,2,1)), |
---|
1587 | & parentarray(1,1,2),parentarray(1,2,2), |
---|
1588 | & parentarray(2,1,2),parentarray(2,2,2), |
---|
1589 | & parentarray(3,1,2),parentarray(3,2,2), |
---|
1590 | & parentarray(4,1,2),parentarray(4,2,2), |
---|
1591 | & parentarray(5,1,2),parentarray(5,2,2), |
---|
1592 | & parentarray(6,1,2),parentarray(6,2,2), |
---|
1593 | & .FALSE.,nbin,ndirin |
---|
1594 | & ) |
---|
1595 | END SELECT |
---|
1596 | CALL Agrif_ParentGrid_to_ChildGrid() |
---|
1597 | ELSE |
---|
1598 | SELECT CASE(nbdim) |
---|
1599 | CASE(1) |
---|
1600 | parent%var%array1(parentarray(1,1,2):parentarray(1,2,2)) = |
---|
1601 | & tempPextend%var%array1( |
---|
1602 | & parentarray(1,1,1):parentarray(1,2,1)) |
---|
1603 | CASE(2) |
---|
1604 | parent%var%array2(parentarray(1,1,2):parentarray(1,2,2), |
---|
1605 | & parentarray(2,1,2):parentarray(2,2,2)) = |
---|
1606 | & tempPextend%var%array2( |
---|
1607 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1608 | & parentarray(2,1,1):parentarray(2,2,1)) |
---|
1609 | CASE(3) |
---|
1610 | parent%var%array3(parentarray(1,1,2):parentarray(1,2,2), |
---|
1611 | & parentarray(2,1,2):parentarray(2,2,2), |
---|
1612 | & parentarray(3,1,2):parentarray(3,2,2)) = |
---|
1613 | & tempPextend%var%array3( |
---|
1614 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1615 | & parentarray(2,1,1):parentarray(2,2,1), |
---|
1616 | & parentarray(3,1,1):parentarray(3,2,1)) |
---|
1617 | CASE(4) |
---|
1618 | parent%var%array4(parentarray(1,1,2):parentarray(1,2,2), |
---|
1619 | & parentarray(2,1,2):parentarray(2,2,2), |
---|
1620 | & parentarray(3,1,2):parentarray(3,2,2), |
---|
1621 | & parentarray(4,1,2):parentarray(4,2,2)) = |
---|
1622 | & tempPextend%var%array4( |
---|
1623 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1624 | & parentarray(2,1,1):parentarray(2,2,1), |
---|
1625 | & parentarray(3,1,1):parentarray(3,2,1), |
---|
1626 | & parentarray(4,1,1):parentarray(4,2,1)) |
---|
1627 | CASE(5) |
---|
1628 | parent%var%array5(parentarray(1,1,2):parentarray(1,2,2), |
---|
1629 | & parentarray(2,1,2):parentarray(2,2,2), |
---|
1630 | & parentarray(3,1,2):parentarray(3,2,2), |
---|
1631 | & parentarray(4,1,2):parentarray(4,2,2), |
---|
1632 | & parentarray(5,1,2):parentarray(5,2,2)) = |
---|
1633 | & tempPextend%var%array5( |
---|
1634 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1635 | & parentarray(2,1,1):parentarray(2,2,1), |
---|
1636 | & parentarray(3,1,1):parentarray(3,2,1), |
---|
1637 | & parentarray(4,1,1):parentarray(4,2,1), |
---|
1638 | & parentarray(5,1,1):parentarray(5,2,1)) |
---|
1639 | CASE(6) |
---|
1640 | parent%var%array6(parentarray(1,1,2):parentarray(1,2,2), |
---|
1641 | & parentarray(2,1,2):parentarray(2,2,2), |
---|
1642 | & parentarray(3,1,2):parentarray(3,2,2), |
---|
1643 | & parentarray(4,1,2):parentarray(4,2,2), |
---|
1644 | & parentarray(5,1,2):parentarray(5,2,2), |
---|
1645 | & parentarray(6,1,2):parentarray(6,2,2)) = |
---|
1646 | & tempPextend%var%array6( |
---|
1647 | & parentarray(1,1,1):parentarray(1,2,1), |
---|
1648 | & parentarray(2,1,1):parentarray(2,2,1), |
---|
1649 | & parentarray(3,1,1):parentarray(3,2,1), |
---|
1650 | & parentarray(4,1,1):parentarray(4,2,1), |
---|
1651 | & parentarray(5,1,1):parentarray(5,2,1), |
---|
1652 | & parentarray(6,1,1):parentarray(6,2,1)) |
---|
1653 | END SELECT |
---|
1654 | ENDIF |
---|
1655 | |
---|
1656 | Call Agrif_nbdim_deallocation(tempPextend%var,nbdim) |
---|
1657 | ENDIF |
---|
1658 | C |
---|
1659 | C |
---|
1660 | C Deallocations |
---|
1661 | |
---|
1662 | IF (memberin) THEN |
---|
1663 | #ifdef AGRIF_MPI |
---|
1664 | Call Agrif_nbdim_deallocation(tempP%var,nbdim) |
---|
1665 | Call Agrif_nbdim_deallocation(tempC%var,nbdim) |
---|
1666 | ! Deallocate(tempC % var) |
---|
1667 | #endif |
---|
1668 | ! Deallocate(tempP % var) |
---|
1669 | ENDIF |
---|
1670 | #ifdef AGRIF_MPI |
---|
1671 | ! Deallocate(tempPextend%var) |
---|
1672 | ! IF (.Not.memberin) Deallocate(tempCextend%var) |
---|
1673 | #endif |
---|
1674 | |
---|
1675 | C |
---|
1676 | C |
---|
1677 | End Subroutine Agrif_UpdatenD |
---|
1678 | C |
---|
1679 | C |
---|
1680 | C ************************************************************************** |
---|
1681 | CCC Subroutine Agrif_Prtbounds |
---|
1682 | C ************************************************************************** |
---|
1683 | C |
---|
1684 | Subroutine Agrif_Prtbounds(nbdim,indmin,indmax,s_Parent_temp, |
---|
1685 | & s_Child_temp,s_Child,ds_Child, |
---|
1686 | & s_Parent,ds_Parent, |
---|
1687 | & pttruetab,cetruetab,pttab_Child, |
---|
1688 | & pttab_Parent, |
---|
1689 | & posvartab_child,TypeUpdate, |
---|
1690 | & loctab_Child |
---|
1691 | #ifdef AGRIF_MPI |
---|
1692 | & ,pttruetabwhole,cetruetabwhole |
---|
1693 | #endif |
---|
1694 | & ) |
---|
1695 | C |
---|
1696 | CCC Description: |
---|
1697 | CCC Subroutine calculating the bounds of the parent grid to be updated |
---|
1698 | CCC by the child grid |
---|
1699 | C |
---|
1700 | C |
---|
1701 | C Declarations: |
---|
1702 | C |
---|
1703 | |
---|
1704 | C |
---|
1705 | #ifdef AGRIF_MPI |
---|
1706 | cccccccccccccccccccccccccc#include "mpif.h" |
---|
1707 | #endif |
---|
1708 | C |
---|
1709 | C Arguments |
---|
1710 | INTEGER :: nbdim |
---|
1711 | INTEGER,DIMENSION(nbdim) :: indmin,indmax |
---|
1712 | REAL,DIMENSION(nbdim) :: s_Parent_temp,s_child_temp |
---|
1713 | REAL,DIMENSION(nbdim) :: s_Child,ds_child |
---|
1714 | REAL,DIMENSION(nbdim) :: s_Parent,ds_Parent |
---|
1715 | INTEGER,DIMENSION(nbdim) :: pttruetab,cetruetab |
---|
1716 | INTEGER,DIMENSION(nbdim) :: posvartab_child,TypeUpdate |
---|
1717 | INTEGER,DIMENSION(nbdim) :: loctab_Child |
---|
1718 | INTEGER,DIMENSION(nbdim) :: pttab_Child,pttab_Parent |
---|
1719 | C |
---|
1720 | C Local variables |
---|
1721 | INTEGER :: i |
---|
1722 | REAL,DIMENSION(nbdim) :: dim_newmin,dim_newmax |
---|
1723 | #ifdef AGRIF_MPI |
---|
1724 | INTEGER,DIMENSION(nbdim) :: pttruetabwhole,cetruetabwhole |
---|
1725 | REAL :: positionmin,positionmax |
---|
1726 | INTEGER :: imin,imax |
---|
1727 | INTEGER :: coeffraf |
---|
1728 | #endif |
---|
1729 | C |
---|
1730 | C |
---|
1731 | do i = 1,nbdim |
---|
1732 | C |
---|
1733 | dim_newmin(i) = s_Child(i) + (pttruetab(i) - |
---|
1734 | & pttab_Child(i)) * ds_Child(i) |
---|
1735 | C |
---|
1736 | dim_newmax(i) = s_Child(i) + (cetruetab(i) - |
---|
1737 | & pttab_Child(i)) * ds_Child(i) |
---|
1738 | C |
---|
1739 | indmin(i) = pttab_Parent(i) + |
---|
1740 | & agrif_ceiling((dim_newmin(i)-s_Parent(i))/ds_Parent(i)) |
---|
1741 | C |
---|
1742 | indmax(i) = pttab_Parent(i) + |
---|
1743 | & agrif_int((dim_newmax(i)-s_Parent(i))/ds_Parent(i)) |
---|
1744 | C |
---|
1745 | #ifdef AGRIF_MPI |
---|
1746 | positionmin = s_Parent(i) + (indmin(i)- |
---|
1747 | & pttab_Parent(i))*ds_Parent(i) |
---|
1748 | IF (loctab_Child(i) .NE. -3) THEN |
---|
1749 | IF (posvartab_child(i) == 1) THEN |
---|
1750 | IF (TypeUpdate(i) .EQ. Agrif_Update_Average) THEN |
---|
1751 | positionmin = positionmin - ds_Parent(i)/2. |
---|
1752 | ELSE IF (TypeUpdate(i) .EQ. Agrif_Update_Full_Weighting) THEN |
---|
1753 | positionmin = positionmin - (ds_Parent(i)-ds_Child(i)) |
---|
1754 | ENDIF |
---|
1755 | ELSE |
---|
1756 | IF (TypeUpdate(i).NE.Agrif_Update_Full_Weighting) THEN |
---|
1757 | positionmin = positionmin - ds_Parent(i)/2. |
---|
1758 | ELSE |
---|
1759 | coeffraf = nint(ds_Parent(i)/ds_Child(i)) |
---|
1760 | if (mod(coeffraf,2) == 1) then |
---|
1761 | positionmin = positionmin - (ds_Parent(i)-ds_Child(i)) |
---|
1762 | else |
---|
1763 | positionmin = positionmin - (ds_Parent(i)-ds_Child(i))-ds_Child(i)/2. |
---|
1764 | endif |
---|
1765 | ENDIF |
---|
1766 | ENDIF |
---|
1767 | ENDIF |
---|
1768 | imin = pttab_Child(i) + |
---|
1769 | & agrif_ceiling((positionmin-s_Child(i))/ds_Child(i)) |
---|
1770 | |
---|
1771 | positionmin = s_Child(i) + (imin - |
---|
1772 | & pttab_Child(i)) * ds_Child(i) |
---|
1773 | |
---|
1774 | pttruetabwhole(i) = imin |
---|
1775 | |
---|
1776 | positionmax = s_Parent(i) + (indmax(i)- |
---|
1777 | & pttab_Parent(i))*ds_Parent(i) |
---|
1778 | IF (loctab_Child(i) .NE. -3) THEN |
---|
1779 | IF (posvartab_child(i) == 1) THEN |
---|
1780 | IF (TypeUpdate(i) .EQ. Agrif_Update_Average) THEN |
---|
1781 | positionmax = positionmax + ds_Parent(i)/2. |
---|
1782 | ELSE IF (TypeUpdate(i) .EQ. Agrif_Update_Full_Weighting) THEN |
---|
1783 | positionmax = positionmax + (ds_Parent(i)-ds_Child(i)) |
---|
1784 | ENDIF |
---|
1785 | ELSE |
---|
1786 | IF (TypeUpdate(i).NE.Agrif_Update_Full_Weighting) THEN |
---|
1787 | positionmax = positionmax + ds_Parent(i)/2. |
---|
1788 | ELSE |
---|
1789 | coeffraf = nint(ds_Parent(i)/ds_Child(i)) |
---|
1790 | if (mod(coeffraf,2) == 1) then |
---|
1791 | positionmax = positionmax + (ds_Parent(i)-ds_Child(i)) |
---|
1792 | else |
---|
1793 | positionmax = positionmax + (ds_Parent(i)-ds_Child(i)) + ds_Child(i)/2. |
---|
1794 | endif |
---|
1795 | |
---|
1796 | ENDIF |
---|
1797 | ENDIF |
---|
1798 | ENDIF |
---|
1799 | imax = pttab_Child(i) + |
---|
1800 | & agrif_int((positionmax-s_Child(i))/ds_Child(i)) |
---|
1801 | |
---|
1802 | positionmax = s_Child(i) + (imax - |
---|
1803 | & pttab_Child(i)) * ds_Child(i) |
---|
1804 | |
---|
1805 | cetruetabwhole(i) = imax |
---|
1806 | |
---|
1807 | #endif |
---|
1808 | C |
---|
1809 | s_Parent_temp(i) = s_Parent(i) + |
---|
1810 | & (indmin(i) - pttab_Parent(i)) * |
---|
1811 | & ds_Parent(i) |
---|
1812 | C |
---|
1813 | s_Child_temp(i) = dim_newmin(i) |
---|
1814 | |
---|
1815 | #ifdef AGRIF_MPI |
---|
1816 | s_Child_temp(i) = positionmin |
---|
1817 | #endif |
---|
1818 | C |
---|
1819 | enddo |
---|
1820 | C |
---|
1821 | Return |
---|
1822 | C |
---|
1823 | C |
---|
1824 | End Subroutine Agrif_Prtbounds |
---|
1825 | C |
---|
1826 | C |
---|
1827 | C |
---|
1828 | C |
---|
1829 | |
---|
1830 | C |
---|
1831 | C |
---|
1832 | C |
---|
1833 | C ************************************************************************** |
---|
1834 | CCC Subroutine Agrif_Update_2D_Recursive |
---|
1835 | C ************************************************************************** |
---|
1836 | C |
---|
1837 | Subroutine Agrif_Update_2D_recursive(TypeUpdate,tempP,tempC, |
---|
1838 | & indmin,indmax, |
---|
1839 | & pttab_child,petab_child, |
---|
1840 | & s_child,s_parent, |
---|
1841 | & ds_child,ds_parent,nbdim) |
---|
1842 | C |
---|
1843 | CCC Description: |
---|
1844 | CCC Subroutine to update a 2D grid variable on the parent grid. |
---|
1845 | CCC It calls Agrif_Update_1D_Recursive and Agrif_UpdateBase. |
---|
1846 | C |
---|
1847 | CC Method: |
---|
1848 | C |
---|
1849 | C Declarations: |
---|
1850 | C |
---|
1851 | |
---|
1852 | C |
---|
1853 | INTEGER :: nbdim |
---|
1854 | INTEGER, DIMENSION(nbdim) :: TypeUpdate ! TYPE of update (copy or average) |
---|
1855 | INTEGER, DIMENSION(nbdim) :: indmin,indmax |
---|
1856 | INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child |
---|
1857 | REAL, DIMENSION(nbdim) :: s_child,s_parent |
---|
1858 | REAL, DIMENSION(nbdim) :: ds_child,ds_parent |
---|
1859 | REAL, DIMENSION(indmin(1):indmax(1), |
---|
1860 | & indmin(2):indmax(2)) :: tempP |
---|
1861 | C REAL, DIMENSION(pttab_child(1):petab_child(1), |
---|
1862 | C & pttab_child(2):petab_child(2)) :: tempC |
---|
1863 | |
---|
1864 | REAL, DIMENSION(:,:) :: tempC |
---|
1865 | C |
---|
1866 | C Local variables |
---|
1867 | REAL, DIMENSION(indmin(1):indmax(1), |
---|
1868 | & pttab_child(2):petab_child(2)) :: tabtemp |
---|
1869 | REAL, DIMENSION(indmin(2):indmax(2), |
---|
1870 | & indmin(1):indmax(1)) :: tempP_trsp |
---|
1871 | REAL, DIMENSION(pttab_child(2):petab_child(2), |
---|
1872 | & indmin(1):indmax(1)) :: tabtemp_trsp |
---|
1873 | INTEGER :: i,j |
---|
1874 | INTEGER :: coeffraf,locind_child_left |
---|
1875 | C |
---|
1876 | tabtemp = 0. |
---|
1877 | |
---|
1878 | |
---|
1879 | coeffraf = nint ( ds_parent(1) / ds_child(1) ) |
---|
1880 | IF((TypeUpdate(1) == AGRIF_Update_average) |
---|
1881 | & .AND. (coeffraf /= 1 ))THEN |
---|
1882 | !---CDIR NEXPAND |
---|
1883 | IF(.NOT. precomputedone(1) ) Call average1Dprecompute2D |
---|
1884 | & (petab_child(2)-pttab_child(2)+1, |
---|
1885 | & indmax(1)-indmin(1)+1, |
---|
1886 | & petab_child(1)-pttab_child(1)+1, |
---|
1887 | & s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
1888 | !---CDIR NEXPAND |
---|
1889 | Call average1Daftercompute |
---|
1890 | & ( tabtemp, tempC, |
---|
1891 | & size(tabtemp), size(tempC), |
---|
1892 | & s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
1893 | |
---|
1894 | ELSEIF((TypeUpdate(1) == AGRIF_Update_copy) |
---|
1895 | & .AND. (coeffraf /= 1 ))THEN |
---|
1896 | !---CDIR NEXPAND |
---|
1897 | IF(.NOT. precomputedone(1) ) Call copy1Dprecompute2D |
---|
1898 | & (petab_child(2)-pttab_child(2)+1, |
---|
1899 | & indmax(1)-indmin(1)+1, |
---|
1900 | & petab_child(1)-pttab_child(1)+1, |
---|
1901 | & s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
1902 | !---CDIR NEXPAND |
---|
1903 | Call copy1Daftercompute |
---|
1904 | & ( tabtemp, tempC, |
---|
1905 | & size(tabtemp), size(tempC), |
---|
1906 | & s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
1907 | |
---|
1908 | ELSE |
---|
1909 | do j = pttab_child(nbdim),petab_child(nbdim) |
---|
1910 | C |
---|
1911 | !---CDIR NEXPAND |
---|
1912 | Call Agrif_Update_1D_recursive(TypeUpdate(1:nbdim-1), |
---|
1913 | & tabtemp(:,j), |
---|
1914 | & tempC(:,j-pttab_child(nbdim)+1), |
---|
1915 | & indmin(1:nbdim-1),indmax(1:nbdim-1), |
---|
1916 | & pttab_child(1:nbdim-1),petab_child(1:nbdim-1), |
---|
1917 | & s_child(1:nbdim-1),s_parent(1:nbdim-1), |
---|
1918 | & ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1) |
---|
1919 | C |
---|
1920 | enddo |
---|
1921 | ENDIF |
---|
1922 | tabtemp_trsp = TRANSPOSE(tabtemp) |
---|
1923 | coeffraf = nint(ds_parent(nbdim)/ds_child(nbdim)) |
---|
1924 | |
---|
1925 | !---CDIR NEXPAND |
---|
1926 | Call Agrif_Compute_nbdim_update(s_parent(nbdim),s_child(nbdim), |
---|
1927 | & ds_parent(nbdim),ds_child(nbdim),coeffraf,locind_child_left) |
---|
1928 | C |
---|
1929 | |
---|
1930 | tempP_trsp = 0. |
---|
1931 | |
---|
1932 | IF((TypeUpdate(2) == AGRIF_Update_average) |
---|
1933 | & .AND. (coeffraf /= 1 ))THEN |
---|
1934 | !---CDIR NEXPAND |
---|
1935 | IF(.NOT. precomputedone(2) ) Call average1Dprecompute2D |
---|
1936 | & ( indmax(1)-indmin(1)+1, |
---|
1937 | & indmax(2)-indmin(2)+1, |
---|
1938 | & petab_child(2)-pttab_child(2)+1, |
---|
1939 | & s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
1940 | !---CDIR NEXPAND |
---|
1941 | Call average1Daftercompute |
---|
1942 | & ( tempP_trsp, tabtemp_trsp, |
---|
1943 | & size(tempP_trsp), size(tabtemp_trsp), |
---|
1944 | & s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
1945 | |
---|
1946 | ELSEIF((TypeUpdate(2) == AGRIF_Update_copy) |
---|
1947 | & .AND. (coeffraf /= 1 ))THEN |
---|
1948 | !---CDIR NEXPAND |
---|
1949 | IF(.NOT. precomputedone(2) ) Call copy1Dprecompute2D |
---|
1950 | & ( indmax(1)-indmin(1)+1, |
---|
1951 | & indmax(2)-indmin(2)+1, |
---|
1952 | & petab_child(2)-pttab_child(2)+1, |
---|
1953 | & s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
1954 | !---CDIR NEXPAND |
---|
1955 | Call copy1Daftercompute |
---|
1956 | & ( tempP_trsp, tabtemp_trsp, |
---|
1957 | & size(tempP_trsp), size(tabtemp_trsp), |
---|
1958 | & s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
1959 | |
---|
1960 | ELSE |
---|
1961 | |
---|
1962 | do i = indmin(1),indmax(1) |
---|
1963 | C |
---|
1964 | !---CDIR NEXPAND |
---|
1965 | Call Agrif_UpdateBase(TypeUpdate(2), |
---|
1966 | & tempP_trsp(indmin(nbdim):indmax(nbdim),i), |
---|
1967 | & tabtemp_trsp(pttab_child(nbdim):petab_child(nbdim),i), |
---|
1968 | & indmin(nbdim),indmax(nbdim), |
---|
1969 | & pttab_child(nbdim),petab_child(nbdim), |
---|
1970 | & s_parent(nbdim),s_child(nbdim), |
---|
1971 | & ds_parent(nbdim),ds_child(nbdim), |
---|
1972 | & coeffraf,locind_child_left) |
---|
1973 | C |
---|
1974 | enddo |
---|
1975 | |
---|
1976 | ENDIF |
---|
1977 | |
---|
1978 | tempP = TRANSPOSE(tempP_trsp) |
---|
1979 | C |
---|
1980 | Return |
---|
1981 | C |
---|
1982 | C |
---|
1983 | End Subroutine Agrif_Update_2D_recursive |
---|
1984 | C |
---|
1985 | C |
---|
1986 | C |
---|
1987 | C ************************************************************************** |
---|
1988 | CCC Subroutine Agrif_Update_3D_Recursive |
---|
1989 | C ************************************************************************** |
---|
1990 | C |
---|
1991 | Subroutine Agrif_Update_3D_recursive(TypeUpdate,tempP,tempC, |
---|
1992 | & indmin,indmax, |
---|
1993 | & pttab_child,petab_child, |
---|
1994 | & s_child,s_parent, |
---|
1995 | & ds_child,ds_parent,nbdim) |
---|
1996 | C |
---|
1997 | CCC Description: |
---|
1998 | CCC Subroutine to update a 3D grid variable on the parent grid. |
---|
1999 | CCC It calls Agrif_Update_2D_Recursive and Agrif_UpdateBase. |
---|
2000 | C |
---|
2001 | CC Method: |
---|
2002 | C |
---|
2003 | C Declarations: |
---|
2004 | C |
---|
2005 | |
---|
2006 | C |
---|
2007 | INTEGER :: nbdim |
---|
2008 | INTEGER, DIMENSION(nbdim) :: TypeUpdate ! TYPE of update (copy or average) |
---|
2009 | INTEGER, DIMENSION(nbdim) :: indmin,indmax |
---|
2010 | INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child |
---|
2011 | REAL, DIMENSION(nbdim) :: s_child,s_parent |
---|
2012 | REAL, DIMENSION(nbdim) :: ds_child,ds_parent |
---|
2013 | REAL, DIMENSION(indmin(1):indmax(1), |
---|
2014 | & indmin(2):indmax(2), |
---|
2015 | & indmin(3):indmax(3)) :: tempP |
---|
2016 | REAL, DIMENSION(pttab_child(1):petab_child(1), |
---|
2017 | & pttab_child(2):petab_child(2), |
---|
2018 | & pttab_child(3):petab_child(3)) :: tempC |
---|
2019 | C |
---|
2020 | C Local variables |
---|
2021 | REAL, DIMENSION(indmin(1):indmax(1), |
---|
2022 | & indmin(2):indmax(2), |
---|
2023 | & pttab_child(3):petab_child(3)) :: tabtemp |
---|
2024 | INTEGER :: i,j,k |
---|
2025 | INTEGER :: coeffraf,locind_child_left |
---|
2026 | INTEGER :: kdeb |
---|
2027 | C |
---|
2028 | C |
---|
2029 | coeffraf = nint ( ds_parent(1) / ds_child(1) ) |
---|
2030 | IF((TypeUpdate(1) == AGRIF_Update_average) |
---|
2031 | & .AND. (coeffraf /= 1 ))THEN |
---|
2032 | !---CDIR NEXPAND |
---|
2033 | Call average1Dprecompute2D |
---|
2034 | & (petab_child(2)-pttab_child(2)+1, |
---|
2035 | & indmax(1)-indmin(1)+1, |
---|
2036 | & petab_child(1)-pttab_child(1)+1, |
---|
2037 | & s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
2038 | precomputedone(1) = .TRUE. |
---|
2039 | ELSEIF((TypeUpdate(1) == AGRIF_Update_copy) |
---|
2040 | & .AND. (coeffraf /= 1 ))THEN |
---|
2041 | !---CDIR NEXPAND |
---|
2042 | Call copy1Dprecompute2D |
---|
2043 | & (petab_child(2)-pttab_child(2)+1, |
---|
2044 | & indmax(1)-indmin(1)+1, |
---|
2045 | & petab_child(1)-pttab_child(1)+1, |
---|
2046 | & s_parent(1),s_child(1),ds_parent(1),ds_child(1),1) |
---|
2047 | precomputedone(1) = .TRUE. |
---|
2048 | ENDIF |
---|
2049 | |
---|
2050 | coeffraf = nint ( ds_parent(2) / ds_child(2) ) |
---|
2051 | IF((TypeUpdate(2) == AGRIF_Update_average) |
---|
2052 | & .AND. (coeffraf /= 1 ))THEN |
---|
2053 | !---CDIR NEXPAND |
---|
2054 | Call average1Dprecompute2D |
---|
2055 | & ( indmax(1)-indmin(1)+1, |
---|
2056 | & indmax(2)-indmin(2)+1, |
---|
2057 | & petab_child(2)-pttab_child(2)+1, |
---|
2058 | & s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
2059 | precomputedone(2) = .TRUE. |
---|
2060 | ELSEIF((TypeUpdate(2) == AGRIF_Update_copy) |
---|
2061 | & .AND. (coeffraf /= 1 ))THEN |
---|
2062 | !---CDIR NEXPAND |
---|
2063 | Call copy1Dprecompute2D |
---|
2064 | & ( indmax(1)-indmin(1)+1, |
---|
2065 | & indmax(2)-indmin(2)+1, |
---|
2066 | & petab_child(2)-pttab_child(2)+1, |
---|
2067 | & s_parent(2),s_child(2),ds_parent(2),ds_child(2),2) |
---|
2068 | precomputedone(2) = .TRUE. |
---|
2069 | ENDIF |
---|
2070 | |
---|
2071 | |
---|
2072 | do k = pttab_child(nbdim),petab_child(nbdim) |
---|
2073 | C |
---|
2074 | Call Agrif_Update_2D_recursive(TypeUpdate(1:nbdim-1), |
---|
2075 | & tabtemp(:,:,k), |
---|
2076 | & tempC(:,:,k), |
---|
2077 | & indmin(1:nbdim-1),indmax(1:nbdim-1), |
---|
2078 | & pttab_child(1:nbdim-1),petab_child(1:nbdim-1), |
---|
2079 | & s_child(1:nbdim-1),s_parent(1:nbdim-1), |
---|
2080 | & ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1) |
---|
2081 | C |
---|
2082 | enddo |
---|
2083 | C |
---|
2084 | precomputedone(1) = .FALSE. |
---|
2085 | precomputedone(2) = .FALSE. |
---|
2086 | coeffraf = nint ( ds_parent(3) / ds_child(3) ) |
---|
2087 | |
---|
2088 | Call Agrif_Compute_nbdim_update(s_parent(nbdim),s_child(nbdim), |
---|
2089 | & ds_parent(nbdim),ds_child(nbdim),coeffraf,locind_child_left) |
---|
2090 | |
---|
2091 | IF (coeffraf == 1) THEN |
---|
2092 | |
---|
2093 | kdeb = pttab_child(3)+locind_child_left-2 |
---|
2094 | do k=indmin(3),indmax(3) |
---|
2095 | kdeb = kdeb + 1 |
---|
2096 | do j = indmin(2),indmax(2) |
---|
2097 | do i = indmin(1),indmax(1) |
---|
2098 | tempP(i,j,k) = tabtemp(i,j,kdeb) |
---|
2099 | enddo |
---|
2100 | enddo |
---|
2101 | enddo |
---|
2102 | |
---|
2103 | ELSE |
---|
2104 | tempP = 0. |
---|
2105 | C |
---|
2106 | do j = indmin(2),indmax(2) |
---|
2107 | C |
---|
2108 | do i = indmin(1),indmax(1) |
---|
2109 | C |
---|
2110 | Call Agrif_UpdateBase(TypeUpdate(3), |
---|
2111 | & tempP(i,j,:), |
---|
2112 | & tabtemp(i,j,:), |
---|
2113 | & indmin(nbdim),indmax(nbdim), |
---|
2114 | & pttab_child(nbdim),petab_child(nbdim), |
---|
2115 | & s_parent(nbdim),s_child(nbdim), |
---|
2116 | & ds_parent(nbdim),ds_child(nbdim), |
---|
2117 | & coeffraf,locind_child_left) |
---|
2118 | C |
---|
2119 | enddo |
---|
2120 | C |
---|
2121 | enddo |
---|
2122 | ENDIF |
---|
2123 | C |
---|
2124 | Return |
---|
2125 | C |
---|
2126 | C |
---|
2127 | End Subroutine Agrif_Update_3D_recursive |
---|
2128 | C |
---|
2129 | C |
---|
2130 | C |
---|
2131 | C ************************************************************************** |
---|
2132 | CCC Subroutine Agrif_Update_4D_Recursive |
---|
2133 | C ************************************************************************** |
---|
2134 | C |
---|
2135 | Subroutine Agrif_Update_4D_recursive(TypeUpdate,tempP,tempC, |
---|
2136 | & indmin,indmax, |
---|
2137 | & pttab_child,petab_child, |
---|
2138 | & s_child,s_parent, |
---|
2139 | & ds_child,ds_parent,nbdim) |
---|
2140 | C |
---|
2141 | CCC Description: |
---|
2142 | CCC Subroutine to update a 4D grid variable on the parent grid. |
---|
2143 | CCC It calls Agrif_Update_3D_Recursive and Agrif_UpdateBase. |
---|
2144 | C |
---|
2145 | CC Method: |
---|
2146 | C |
---|
2147 | C Declarations: |
---|
2148 | C |
---|
2149 | |
---|
2150 | C |
---|
2151 | INTEGER :: nbdim |
---|
2152 | INTEGER, DIMENSION(nbdim) :: TypeUpdate ! TYPE of update (copy or average) |
---|
2153 | INTEGER, DIMENSION(nbdim) :: indmin,indmax |
---|
2154 | INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child |
---|
2155 | REAL, DIMENSION(nbdim) :: s_child,s_parent |
---|
2156 | REAL, DIMENSION(nbdim) :: ds_child,ds_parent |
---|
2157 | REAL, DIMENSION(indmin(1):indmax(1), |
---|
2158 | & indmin(2):indmax(2), |
---|
2159 | & indmin(3):indmax(3), |
---|
2160 | & indmin(4):indmax(4)) :: tempP |
---|
2161 | REAL, DIMENSION(pttab_child(1):petab_child(1), |
---|
2162 | & pttab_child(2):petab_child(2), |
---|
2163 | & pttab_child(3):petab_child(3), |
---|
2164 | & pttab_child(4):petab_child(4)) :: tempC |
---|
2165 | C |
---|
2166 | C Local variables |
---|
2167 | REAL, DIMENSION(:,:,:,:), Allocatable :: tabtemp |
---|
2168 | INTEGER :: i,j,k,l |
---|
2169 | INTEGER :: coeffraf,locind_child_left |
---|
2170 | C |
---|
2171 | C |
---|
2172 | Allocate(tabtemp(indmin(1):indmax(1), |
---|
2173 | & indmin(2):indmax(2), |
---|
2174 | & indmin(3):indmax(3), |
---|
2175 | & pttab_child(4):petab_child(4))) |
---|
2176 | C |
---|
2177 | do l = pttab_child(nbdim),petab_child(nbdim) |
---|
2178 | C |
---|
2179 | Call Agrif_Update_3D_recursive(TypeUpdate(1:nbdim-1), |
---|
2180 | & tabtemp(indmin(nbdim-3):indmax(nbdim-3), |
---|
2181 | & indmin(nbdim-2):indmax(nbdim-2), |
---|
2182 | & indmin(nbdim-1):indmax(nbdim-1),l), |
---|
2183 | & tempC(pttab_child(nbdim-3):petab_child(nbdim-3), |
---|
2184 | & pttab_child(nbdim-2):petab_child(nbdim-2), |
---|
2185 | & pttab_child(nbdim-1):petab_child(nbdim-1),l), |
---|
2186 | & indmin(1:nbdim-1),indmax(1:nbdim-1), |
---|
2187 | & pttab_child(1:nbdim-1),petab_child(1:nbdim-1), |
---|
2188 | & s_child(1:nbdim-1),s_parent(1:nbdim-1), |
---|
2189 | & ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1) |
---|
2190 | C |
---|
2191 | enddo |
---|
2192 | |
---|
2193 | Call Agrif_Compute_nbdim_update(s_parent(nbdim),s_child(nbdim), |
---|
2194 | & ds_parent(nbdim),ds_child(nbdim),coeffraf,locind_child_left) |
---|
2195 | C |
---|
2196 | tempP = 0. |
---|
2197 | |
---|
2198 | do k = indmin(3),indmax(3) |
---|
2199 | C |
---|
2200 | do j = indmin(2),indmax(2) |
---|
2201 | C |
---|
2202 | do i = indmin(1),indmax(1) |
---|
2203 | C |
---|
2204 | Call Agrif_UpdateBase(TypeUpdate(4), |
---|
2205 | & tempP(i,j,k,indmin(nbdim):indmax(nbdim)), |
---|
2206 | & tabtemp(i,j,k,pttab_child(nbdim):petab_child(nbdim)), |
---|
2207 | & indmin(nbdim),indmax(nbdim), |
---|
2208 | & pttab_child(nbdim),petab_child(nbdim), |
---|
2209 | & s_parent(nbdim),s_child(nbdim), |
---|
2210 | & ds_parent(nbdim),ds_child(nbdim), |
---|
2211 | & coeffraf,locind_child_left) |
---|
2212 | C |
---|
2213 | enddo |
---|
2214 | C |
---|
2215 | enddo |
---|
2216 | C |
---|
2217 | enddo |
---|
2218 | C |
---|
2219 | Deallocate(tabtemp) |
---|
2220 | C |
---|
2221 | Return |
---|
2222 | C |
---|
2223 | C |
---|
2224 | End Subroutine Agrif_Update_4D_recursive |
---|
2225 | C |
---|
2226 | C |
---|
2227 | C |
---|
2228 | C ************************************************************************** |
---|
2229 | CCC Subroutine Agrif_Update_5D_Recursive |
---|
2230 | C ************************************************************************** |
---|
2231 | C |
---|
2232 | Subroutine Agrif_Update_5D_recursive(TypeUpdate,tempP,tempC, |
---|
2233 | & indmin,indmax, |
---|
2234 | & pttab_child,petab_child, |
---|
2235 | & s_child,s_parent, |
---|
2236 | & ds_child,ds_parent,nbdim) |
---|
2237 | C |
---|
2238 | CCC Description: |
---|
2239 | CCC Subroutine to update a 5D grid variable on the parent grid. |
---|
2240 | CCC It calls Agrif_Update_4D_Recursive and Agrif_UpdateBase. |
---|
2241 | C |
---|
2242 | CC Method: |
---|
2243 | C |
---|
2244 | C Declarations: |
---|
2245 | C |
---|
2246 | |
---|
2247 | C |
---|
2248 | INTEGER :: nbdim |
---|
2249 | INTEGER, DIMENSION(nbdim) :: TypeUpdate ! TYPE of update (copy or average) |
---|
2250 | INTEGER, DIMENSION(nbdim) :: indmin,indmax |
---|
2251 | INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child |
---|
2252 | REAL, DIMENSION(nbdim) :: s_child,s_parent |
---|
2253 | REAL, DIMENSION(nbdim) :: ds_child,ds_parent |
---|
2254 | REAL, DIMENSION(indmin(1):indmax(1), |
---|
2255 | & indmin(2):indmax(2), |
---|
2256 | & indmin(3):indmax(3), |
---|
2257 | & indmin(4):indmax(4), |
---|
2258 | & indmin(5):indmax(5)) :: tempP |
---|
2259 | REAL, DIMENSION(pttab_child(1):petab_child(1), |
---|
2260 | & pttab_child(2):petab_child(2), |
---|
2261 | & pttab_child(3):petab_child(3), |
---|
2262 | & pttab_child(4):petab_child(4), |
---|
2263 | & pttab_child(5):petab_child(5)) :: tempC |
---|
2264 | C |
---|
2265 | C Local variables |
---|
2266 | REAL, DIMENSION(:,:,:,:,:), Allocatable :: tabtemp |
---|
2267 | INTEGER :: i,j,k,l,m |
---|
2268 | INTEGER :: coeffraf,locind_child_left |
---|
2269 | C |
---|
2270 | C |
---|
2271 | Allocate(tabtemp(indmin(1):indmax(1), |
---|
2272 | & indmin(2):indmax(2), |
---|
2273 | & indmin(3):indmax(3), |
---|
2274 | & indmin(4):indmax(4), |
---|
2275 | & pttab_child(5):petab_child(5))) |
---|
2276 | C |
---|
2277 | do m = pttab_child(nbdim),petab_child(nbdim) |
---|
2278 | C |
---|
2279 | Call Agrif_Update_4D_recursive(TypeUpdate(1:nbdim-1), |
---|
2280 | & tabtemp(indmin(nbdim-4):indmax(nbdim-4), |
---|
2281 | & indmin(nbdim-3):indmax(nbdim-3), |
---|
2282 | & indmin(nbdim-2):indmax(nbdim-2), |
---|
2283 | & indmin(nbdim-1):indmax(nbdim-1),m), |
---|
2284 | & tempC(pttab_child(nbdim-4):petab_child(nbdim-4), |
---|
2285 | & pttab_child(nbdim-3):petab_child(nbdim-3), |
---|
2286 | & pttab_child(nbdim-2):petab_child(nbdim-2), |
---|
2287 | & pttab_child(nbdim-1):petab_child(nbdim-1),m), |
---|
2288 | & indmin(1:nbdim-1),indmax(1:nbdim-1), |
---|
2289 | & pttab_child(1:nbdim-1),petab_child(1:nbdim-1), |
---|
2290 | & s_child(1:nbdim-1),s_parent(1:nbdim-1), |
---|
2291 | & ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1) |
---|
2292 | C |
---|
2293 | enddo |
---|
2294 | |
---|
2295 | Call Agrif_Compute_nbdim_update(s_parent(nbdim),s_child(nbdim), |
---|
2296 | & ds_parent(nbdim),ds_child(nbdim),coeffraf,locind_child_left) |
---|
2297 | tempP = 0. |
---|
2298 | C |
---|
2299 | do l = indmin(4),indmax(4) |
---|
2300 | C |
---|
2301 | do k = indmin(3),indmax(3) |
---|
2302 | C |
---|
2303 | do j = indmin(2),indmax(2) |
---|
2304 | C |
---|
2305 | do i = indmin(1),indmax(1) |
---|
2306 | C |
---|
2307 | Call Agrif_UpdateBase(TypeUpdate(5), |
---|
2308 | & tempP(i,j,k,l,indmin(nbdim):indmax(nbdim)), |
---|
2309 | & tabtemp(i,j,k,l, |
---|
2310 | & pttab_child(nbdim):petab_child(nbdim)), |
---|
2311 | & indmin(nbdim),indmax(nbdim), |
---|
2312 | & pttab_child(nbdim),petab_child(nbdim), |
---|
2313 | & s_parent(nbdim),s_child(nbdim), |
---|
2314 | & ds_parent(nbdim),ds_child(nbdim), |
---|
2315 | & coeffraf,locind_child_left) |
---|
2316 | C |
---|
2317 | enddo |
---|
2318 | C |
---|
2319 | enddo |
---|
2320 | C |
---|
2321 | enddo |
---|
2322 | C |
---|
2323 | enddo |
---|
2324 | C |
---|
2325 | Deallocate(tabtemp) |
---|
2326 | C |
---|
2327 | Return |
---|
2328 | C |
---|
2329 | C |
---|
2330 | End Subroutine Agrif_Update_5D_recursive |
---|
2331 | C |
---|
2332 | C |
---|
2333 | C |
---|
2334 | C |
---|
2335 | C ************************************************************************** |
---|
2336 | CCC Subroutine Agrif_Update_6D_Recursive |
---|
2337 | C ************************************************************************** |
---|
2338 | C |
---|
2339 | Subroutine Agrif_Update_6D_recursive(TypeUpdate,tempP,tempC, |
---|
2340 | & indmin,indmax, |
---|
2341 | & pttab_child,petab_child, |
---|
2342 | & s_child,s_parent, |
---|
2343 | & ds_child,ds_parent,nbdim) |
---|
2344 | C |
---|
2345 | CCC Description: |
---|
2346 | CCC Subroutine to update a 6D grid variable on the parent grid. |
---|
2347 | CCC It calls Agrif_Update_5D_Recursive and Agrif_UpdateBase. |
---|
2348 | C |
---|
2349 | CC Method: |
---|
2350 | C |
---|
2351 | C Declarations: |
---|
2352 | C |
---|
2353 | |
---|
2354 | C |
---|
2355 | INTEGER :: nbdim |
---|
2356 | INTEGER, DIMENSION(nbdim) :: TypeUpdate ! TYPE of update (copy or average) |
---|
2357 | INTEGER, DIMENSION(nbdim) :: indmin,indmax |
---|
2358 | INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child |
---|
2359 | REAL, DIMENSION(nbdim) :: s_child,s_parent |
---|
2360 | REAL, DIMENSION(nbdim) :: ds_child,ds_parent |
---|
2361 | REAL, DIMENSION(indmin(1):indmax(1), |
---|
2362 | & indmin(2):indmax(2), |
---|
2363 | & indmin(3):indmax(3), |
---|
2364 | & indmin(4):indmax(4), |
---|
2365 | & indmin(5):indmax(5), |
---|
2366 | & indmin(6):indmax(6)) :: tempP |
---|
2367 | REAL, DIMENSION(pttab_child(1):petab_child(1), |
---|
2368 | & pttab_child(2):petab_child(2), |
---|
2369 | & pttab_child(3):petab_child(3), |
---|
2370 | & pttab_child(4):petab_child(4), |
---|
2371 | & pttab_child(5):petab_child(5), |
---|
2372 | & pttab_child(6):petab_child(6)) :: tempC |
---|
2373 | C |
---|
2374 | C Local variables |
---|
2375 | REAL, DIMENSION(:,:,:,:,:,:), Allocatable :: tabtemp |
---|
2376 | INTEGER :: i,j,k,l,m,n |
---|
2377 | INTEGER :: coeffraf,locind_child_left |
---|
2378 | C |
---|
2379 | C |
---|
2380 | Allocate(tabtemp(indmin(1):indmax(1), |
---|
2381 | & indmin(2):indmax(2), |
---|
2382 | & indmin(3):indmax(3), |
---|
2383 | & indmin(4):indmax(4), |
---|
2384 | & indmin(5):indmax(5), |
---|
2385 | & pttab_child(6):petab_child(6))) |
---|
2386 | C |
---|
2387 | do n = pttab_child(nbdim),petab_child(nbdim) |
---|
2388 | C |
---|
2389 | Call Agrif_Update_5D_recursive(TypeUpdate(1:nbdim-1), |
---|
2390 | & tabtemp(indmin(nbdim-5):indmax(nbdim-5), |
---|
2391 | & indmin(nbdim-4):indmax(nbdim-4), |
---|
2392 | & indmin(nbdim-3):indmax(nbdim-3), |
---|
2393 | & indmin(nbdim-2):indmax(nbdim-2), |
---|
2394 | & indmin(nbdim-1):indmax(nbdim-1),n), |
---|
2395 | & tempC(pttab_child(nbdim-5):petab_child(nbdim-5), |
---|
2396 | & pttab_child(nbdim-4):petab_child(nbdim-4), |
---|
2397 | & pttab_child(nbdim-3):petab_child(nbdim-3), |
---|
2398 | & pttab_child(nbdim-2):petab_child(nbdim-2), |
---|
2399 | & pttab_child(nbdim-1):petab_child(nbdim-1),n), |
---|
2400 | & indmin(1:nbdim-1),indmax(1:nbdim-1), |
---|
2401 | & pttab_child(1:nbdim-1),petab_child(1:nbdim-1), |
---|
2402 | & s_child(1:nbdim-1),s_parent(1:nbdim-1), |
---|
2403 | & ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1) |
---|
2404 | C |
---|
2405 | enddo |
---|
2406 | |
---|
2407 | Call Agrif_Compute_nbdim_update(s_parent(nbdim),s_child(nbdim), |
---|
2408 | & ds_parent(nbdim),ds_child(nbdim),coeffraf,locind_child_left) |
---|
2409 | C |
---|
2410 | tempP = 0. |
---|
2411 | |
---|
2412 | do m = indmin(5),indmax(5) |
---|
2413 | do l = indmin(4),indmax(4) |
---|
2414 | C |
---|
2415 | do k = indmin(3),indmax(3) |
---|
2416 | C |
---|
2417 | do j = indmin(2),indmax(2) |
---|
2418 | C |
---|
2419 | do i = indmin(1),indmax(1) |
---|
2420 | C |
---|
2421 | Call Agrif_UpdateBase(TypeUpdate(6), |
---|
2422 | & tempP(i,j,k,l,m,indmin(nbdim):indmax(nbdim)), |
---|
2423 | & tabtemp(i,j,k,l,m, |
---|
2424 | & pttab_child(nbdim):petab_child(nbdim)), |
---|
2425 | & indmin(nbdim),indmax(nbdim), |
---|
2426 | & pttab_child(nbdim),petab_child(nbdim), |
---|
2427 | & s_parent(nbdim),s_child(nbdim), |
---|
2428 | & ds_parent(nbdim),ds_child(nbdim), |
---|
2429 | & coeffraf,locind_child_left) |
---|
2430 | C |
---|
2431 | enddo |
---|
2432 | C |
---|
2433 | enddo |
---|
2434 | C |
---|
2435 | enddo |
---|
2436 | C |
---|
2437 | enddo |
---|
2438 | enddo |
---|
2439 | C |
---|
2440 | Deallocate(tabtemp) |
---|
2441 | C |
---|
2442 | Return |
---|
2443 | C |
---|
2444 | C |
---|
2445 | End Subroutine Agrif_Update_6D_recursive |
---|
2446 | C |
---|
2447 | C |
---|
2448 | C |
---|
2449 | C ************************************************************************** |
---|
2450 | CCC Subroutine Agrif_UpdateBase |
---|
2451 | C ************************************************************************** |
---|
2452 | C |
---|
2453 | Subroutine Agrif_UpdateBase(TypeUpdate, |
---|
2454 | & parenttab,childtab, |
---|
2455 | & indmin,indmax,pttab_child,petab_child, |
---|
2456 | & s_parent,s_child,ds_parent,ds_child, |
---|
2457 | & coeffraf,locind_child_left) |
---|
2458 | C |
---|
2459 | CCC Description: |
---|
2460 | CCC Subroutine calling the updating method chosen by the user (copy, average |
---|
2461 | CCC or full-weighting). |
---|
2462 | C |
---|
2463 | CC Method: |
---|
2464 | C |
---|
2465 | C Declarations: |
---|
2466 | C |
---|
2467 | |
---|
2468 | C |
---|
2469 | INTEGER :: TypeUpdate |
---|
2470 | INTEGER :: indmin,indmax |
---|
2471 | INTEGER :: pttab_child,petab_child |
---|
2472 | REAL,DIMENSION(indmin:indmax) :: parenttab |
---|
2473 | REAL,DIMENSION(pttab_child:petab_child) :: childtab |
---|
2474 | REAL :: s_parent,s_child |
---|
2475 | REAL :: ds_parent,ds_child |
---|
2476 | INTEGER :: coeffraf,locind_child_left |
---|
2477 | C |
---|
2478 | C |
---|
2479 | if (TypeUpdate == AGRIF_Update_copy) then |
---|
2480 | C |
---|
2481 | Call agrif_copy1D |
---|
2482 | & (parenttab,childtab, |
---|
2483 | & indmax-indmin+1,petab_child-pttab_child+1, |
---|
2484 | & s_parent,s_child,ds_parent,ds_child) |
---|
2485 | C |
---|
2486 | elseif (TypeUpdate == AGRIF_Update_average) then |
---|
2487 | C |
---|
2488 | Call average1D |
---|
2489 | & (parenttab,childtab, |
---|
2490 | & indmax-indmin+1,petab_child-pttab_child+1, |
---|
2491 | & s_parent,s_child,ds_parent,ds_child) |
---|
2492 | C |
---|
2493 | elseif (TypeUpdate == AGRIF_Update_full_weighting) then |
---|
2494 | C |
---|
2495 | Call full_weighting1D |
---|
2496 | & (parenttab,childtab, |
---|
2497 | & indmax-indmin+1,petab_child-pttab_child+1, |
---|
2498 | & s_parent,s_child,ds_parent,ds_child, |
---|
2499 | & coeffraf,locind_child_left) |
---|
2500 | C |
---|
2501 | endif |
---|
2502 | C |
---|
2503 | Return |
---|
2504 | C |
---|
2505 | C |
---|
2506 | End Subroutine Agrif_UpdateBase |
---|
2507 | C |
---|
2508 | C |
---|
2509 | |
---|
2510 | Subroutine Agrif_Compute_nbdim_update(s_parent,s_child, |
---|
2511 | & ds_parent,ds_child,coeffraf,locind_child_left) |
---|
2512 | real :: s_parent,s_child,ds_parent,ds_child |
---|
2513 | integer :: coeffraf,locind_child_left |
---|
2514 | |
---|
2515 | coeffraf = nint(ds_parent/ds_child) |
---|
2516 | locind_child_left = 1 + agrif_int((s_parent-s_child)/ds_child) |
---|
2517 | |
---|
2518 | End Subroutine Agrif_Compute_nbdim_update |
---|
2519 | |
---|
2520 | #if defined AGRIF_MPI |
---|
2521 | Subroutine Agrif_Find_list_update(list_update,pttab,petab, |
---|
2522 | & pttab_Child,pttab_Parent,nbdim, |
---|
2523 | & find_list_update,tab4t,tab5t,memberinall,memberinall2, |
---|
2524 | & sendtoproc1,recvfromproc1,sendtoproc2,recvfromproc2) |
---|
2525 | TYPE(Agrif_List_Interp_Loc), Pointer :: list_update |
---|
2526 | INTEGER :: nbdim |
---|
2527 | INTEGER,DIMENSION(nbdim) :: pttab,petab,pttab_Child,pttab_Parent |
---|
2528 | LOGICAL :: find_list_update |
---|
2529 | Type(Agrif_List_Interp_loc), Pointer :: parcours |
---|
2530 | INTEGER :: i |
---|
2531 | C |
---|
2532 | INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,8) :: tab4t |
---|
2533 | INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,8) :: tab5t |
---|
2534 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: memberinall,memberinall2 |
---|
2535 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: sendtoproc1,recvfromproc1 |
---|
2536 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: sendtoproc2,recvfromproc2 |
---|
2537 | |
---|
2538 | find_list_update = .FALSE. |
---|
2539 | |
---|
2540 | parcours => list_update |
---|
2541 | |
---|
2542 | Find_loop : Do While (associated(parcours)) |
---|
2543 | Do i=1,nbdim |
---|
2544 | IF ((pttab(i) /= parcours%interp_loc%pttab(i)).OR. |
---|
2545 | & (petab(i) /= parcours%interp_loc%petab(i)).OR. |
---|
2546 | & (pttab_child(i) /= parcours%interp_loc%pttab_child(i)).OR. |
---|
2547 | & (pttab_parent(i) /= parcours%interp_loc%pttab_parent(i))) |
---|
2548 | & THEN |
---|
2549 | parcours=>parcours%suiv |
---|
2550 | Cycle Find_loop |
---|
2551 | ENDIF |
---|
2552 | EndDo |
---|
2553 | |
---|
2554 | tab4t = parcours%interp_loc%tab4t(1:nbdim,0:Agrif_Nbprocs-1,1:8) |
---|
2555 | memberinall = parcours%interp_loc%memberinall(0:Agrif_Nbprocs-1) |
---|
2556 | |
---|
2557 | tab5t = parcours%interp_loc%tab5t(1:nbdim,0:Agrif_Nbprocs-1,1:8) |
---|
2558 | memberinall2 = |
---|
2559 | & parcours%interp_loc%memberinall2(0:Agrif_Nbprocs-1) |
---|
2560 | sendtoproc1 = |
---|
2561 | & parcours%interp_loc%sendtoproc1(0:Agrif_Nbprocs-1) |
---|
2562 | recvfromproc1 = |
---|
2563 | & parcours%interp_loc%recvfromproc1(0:Agrif_Nbprocs-1) |
---|
2564 | sendtoproc2 = |
---|
2565 | & parcours%interp_loc%sendtoproc2(0:Agrif_Nbprocs-1) |
---|
2566 | recvfromproc2 = |
---|
2567 | & parcours%interp_loc%recvfromproc2(0:Agrif_Nbprocs-1) |
---|
2568 | |
---|
2569 | find_list_update = .TRUE. |
---|
2570 | Exit Find_loop |
---|
2571 | End Do Find_loop |
---|
2572 | |
---|
2573 | End Subroutine Agrif_Find_list_update |
---|
2574 | |
---|
2575 | Subroutine Agrif_AddTo_list_update(list_update,pttab,petab, |
---|
2576 | & pttab_Child,pttab_Parent,nbdim |
---|
2577 | & ,tab4t,tab5t,memberinall,memberinall2, |
---|
2578 | & sendtoproc1,recvfromproc1,sendtoproc2,recvfromproc2) |
---|
2579 | |
---|
2580 | TYPE(Agrif_List_Interp_Loc), Pointer :: list_update |
---|
2581 | INTEGER :: nbdim |
---|
2582 | INTEGER,DIMENSION(nbdim) :: pttab,petab,pttab_Child,pttab_Parent |
---|
2583 | INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,8) :: tab4t |
---|
2584 | INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,8) :: tab5t |
---|
2585 | LOGICAL,DIMENSION(0:Agrif_Nbprocs-1) :: memberinall, memberinall2 |
---|
2586 | LOGICAL,DIMENSION(0:Agrif_Nbprocs-1) :: sendtoproc1, recvfromproc1 |
---|
2587 | LOGICAL,DIMENSION(0:Agrif_Nbprocs-1) :: sendtoproc2, recvfromproc2 |
---|
2588 | |
---|
2589 | Type(Agrif_List_Interp_loc), Pointer :: parcours |
---|
2590 | |
---|
2591 | Allocate(parcours) |
---|
2592 | Allocate(parcours%interp_loc) |
---|
2593 | |
---|
2594 | parcours%interp_loc%pttab(1:nbdim) = pttab(1:nbdim) |
---|
2595 | parcours%interp_loc%petab(1:nbdim) = petab(1:nbdim) |
---|
2596 | parcours%interp_loc%pttab_child(1:nbdim) = pttab_child(1:nbdim) |
---|
2597 | parcours%interp_loc%pttab_parent(1:nbdim) = pttab_parent(1:nbdim) |
---|
2598 | Allocate(parcours%interp_loc%tab4t(nbdim,0:Agrif_Nbprocs-1,8)) |
---|
2599 | Allocate(parcours%interp_loc%memberinall(0:Agrif_Nbprocs-1)) |
---|
2600 | |
---|
2601 | Allocate(parcours%interp_loc%tab5t(nbdim,0:Agrif_Nbprocs-1,8)) |
---|
2602 | Allocate(parcours%interp_loc%memberinall2(0:Agrif_Nbprocs-1)) |
---|
2603 | Allocate(parcours%interp_loc%sendtoproc1(0:Agrif_Nbprocs-1)) |
---|
2604 | Allocate(parcours%interp_loc%recvfromproc1(0:Agrif_Nbprocs-1)) |
---|
2605 | Allocate(parcours%interp_loc%sendtoproc2(0:Agrif_Nbprocs-1)) |
---|
2606 | Allocate(parcours%interp_loc%recvfromproc2(0:Agrif_Nbprocs-1)) |
---|
2607 | |
---|
2608 | parcours%interp_loc%tab4t=tab4t |
---|
2609 | parcours%interp_loc%memberinall=memberinall |
---|
2610 | |
---|
2611 | parcours%interp_loc%tab5t=tab5t |
---|
2612 | parcours%interp_loc%memberinall2=memberinall2 |
---|
2613 | parcours%interp_loc%sendtoproc1=sendtoproc1 |
---|
2614 | parcours%interp_loc%recvfromproc1=recvfromproc1 |
---|
2615 | parcours%interp_loc%sendtoproc2=sendtoproc2 |
---|
2616 | parcours%interp_loc%recvfromproc2=recvfromproc2 |
---|
2617 | |
---|
2618 | parcours%suiv => list_update |
---|
2619 | |
---|
2620 | list_update => parcours |
---|
2621 | |
---|
2622 | End Subroutine Agrif_Addto_list_update |
---|
2623 | #endif |
---|
2624 | |
---|
2625 | |
---|
2626 | C ************************************************************************** |
---|
2627 | CCC Subroutine Agrif_Update_1D_Recursive |
---|
2628 | C ************************************************************************** |
---|
2629 | C |
---|
2630 | Subroutine Agrif_Update_1D_recursive(TypeUpdate,tempP,tempC, |
---|
2631 | & indmin,indmax, |
---|
2632 | & pttab_child,petab_child, |
---|
2633 | & s_child,s_parent, |
---|
2634 | & ds_child,ds_parent,nbdim) |
---|
2635 | C |
---|
2636 | CCC Description: |
---|
2637 | CCC Subroutine to update a 1D grid variable on the parent grid. |
---|
2638 | C |
---|
2639 | CC Method: |
---|
2640 | C |
---|
2641 | C Declarations: |
---|
2642 | C |
---|
2643 | |
---|
2644 | C |
---|
2645 | C Arguments |
---|
2646 | INTEGER :: nbdim |
---|
2647 | INTEGER, DIMENSION(nbdim) :: TypeUpdate ! TYPE of update (copy or average) |
---|
2648 | INTEGER, DIMENSION(nbdim) :: indmin,indmax |
---|
2649 | INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child |
---|
2650 | REAL, DIMENSION(nbdim) :: s_child,s_parent |
---|
2651 | REAL, DIMENSION(nbdim) :: ds_child,ds_parent |
---|
2652 | REAL, DIMENSION(indmin(nbdim):indmax(nbdim)) :: tempP |
---|
2653 | REAL, DIMENSION(pttab_child(nbdim):petab_child(nbdim)) :: tempC |
---|
2654 | INTEGER :: coeffraf,locind_child_left |
---|
2655 | C |
---|
2656 | C |
---|
2657 | Call Agrif_Compute_nbdim_update(s_parent(nbdim),s_child(nbdim), |
---|
2658 | & ds_parent(nbdim),ds_child(nbdim),coeffraf,locind_child_left) |
---|
2659 | |
---|
2660 | Call Agrif_UpdateBase(TypeUpdate(1), |
---|
2661 | & tempP(indmin(nbdim):indmax(nbdim)), |
---|
2662 | & tempC(pttab_child(nbdim):petab_child(nbdim)), |
---|
2663 | & indmin(nbdim),indmax(nbdim), |
---|
2664 | & pttab_child(nbdim),petab_child(nbdim), |
---|
2665 | & s_parent(nbdim),s_child(nbdim), |
---|
2666 | & ds_parent(nbdim),ds_child(nbdim), |
---|
2667 | & coeffraf,locind_child_left) |
---|
2668 | |
---|
2669 | C |
---|
2670 | Return |
---|
2671 | C |
---|
2672 | C |
---|
2673 | End Subroutine Agrif_Update_1D_recursive |
---|
2674 | |
---|
2675 | End Module Agrif_Update |
---|