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