[1901] | 1 | ! |
---|
| 2 | ! $Id: modbcfunction.F 1200 2008-09-24 13:05:20Z rblod $ |
---|
| 3 | ! |
---|
| 4 | C AGRIF (Adaptive Grid Refinement In Fortran) |
---|
| 5 | C |
---|
| 6 | C Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr) |
---|
| 7 | C Christophe Vouland (Christophe.Vouland@imag.fr) |
---|
| 8 | C |
---|
| 9 | C This program is free software; you can redistribute it and/or modify |
---|
| 10 | C it under the terms of the GNU General Public License as published by |
---|
| 11 | C the Free Software Foundation; either version 2 of the License, or |
---|
| 12 | C (at your option) any later version. |
---|
| 13 | C |
---|
| 14 | C This program is distributed in the hope that it will be useful, |
---|
| 15 | C but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
| 16 | C MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
| 17 | C GNU General Public License for more details. |
---|
| 18 | C |
---|
| 19 | C You should have received a copy of the GNU General Public License |
---|
| 20 | C along with this program; if not, write to the Free Software |
---|
| 21 | C Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
---|
| 22 | C |
---|
| 23 | C |
---|
| 24 | C |
---|
| 25 | CCC Module AGRIF_bcfunction |
---|
| 26 | C |
---|
| 27 | C |
---|
| 28 | Module Agrif_bcfunction |
---|
| 29 | CCC Description: |
---|
| 30 | CCC |
---|
| 31 | C |
---|
| 32 | C Modules used: |
---|
| 33 | C |
---|
| 34 | Use Agrif_Boundary |
---|
| 35 | Use Agrif_Update |
---|
| 36 | Use Agrif_fluxmod |
---|
| 37 | C |
---|
| 38 | IMPLICIT NONE |
---|
| 39 | C |
---|
| 40 | interface Agrif_Bc_variable |
---|
| 41 | module procedure Agrif_Bc_variable0d, |
---|
| 42 | & Agrif_Bc_variable1d, |
---|
| 43 | & Agrif_Bc_variable2d, |
---|
| 44 | & Agrif_Bc_variable3d, |
---|
| 45 | & Agrif_Bc_variable4d, |
---|
| 46 | & Agrif_Bc_variable5d |
---|
| 47 | end interface |
---|
| 48 | C |
---|
| 49 | interface Agrif_Set_Parent |
---|
| 50 | module procedure Agrif_Set_Parent_int, |
---|
| 51 | & Agrif_Set_Parent_real |
---|
| 52 | end interface |
---|
| 53 | C |
---|
| 54 | interface Agrif_Interp_variable |
---|
| 55 | module procedure Agrif_Interp_var0d, |
---|
| 56 | & Agrif_Interp_var1d, |
---|
| 57 | & Agrif_Interp_var2d, |
---|
| 58 | & Agrif_Interp_var3d, |
---|
| 59 | & Agrif_Interp_var4d, |
---|
| 60 | & Agrif_Interp_var5d |
---|
| 61 | end interface |
---|
| 62 | C |
---|
| 63 | interface Agrif_Init_variable |
---|
| 64 | module procedure Agrif_Init_variable0d, |
---|
| 65 | & Agrif_Init_variable1d, |
---|
| 66 | & Agrif_Init_variable2d, |
---|
| 67 | & Agrif_Init_variable3d |
---|
| 68 | end interface |
---|
| 69 | C |
---|
| 70 | interface Agrif_update_variable |
---|
| 71 | module procedure Agrif_update_var0d, |
---|
| 72 | & Agrif_update_var1d, |
---|
| 73 | & Agrif_update_var2d, |
---|
| 74 | & Agrif_update_var3d, |
---|
| 75 | & Agrif_update_var4d, |
---|
| 76 | & Agrif_update_var5d |
---|
| 77 | end interface |
---|
| 78 | C |
---|
| 79 | Contains |
---|
| 80 | C |
---|
| 81 | C ************************************************************************** |
---|
| 82 | CCC Subroutine Agrif_Set_type |
---|
| 83 | C ************************************************************************** |
---|
| 84 | C |
---|
| 85 | Subroutine Agrif_Set_type(tabvarsindic,posvar,point) |
---|
| 86 | C |
---|
| 87 | CCC Description: |
---|
| 88 | CCC To set the TYPE of the variable. |
---|
| 89 | C |
---|
| 90 | C Modules used: |
---|
| 91 | C |
---|
| 92 | |
---|
| 93 | C |
---|
| 94 | C Declarations: |
---|
| 95 | C |
---|
| 96 | C |
---|
| 97 | C |
---|
| 98 | C Arguments |
---|
| 99 | C |
---|
| 100 | INTEGER, DIMENSION(:) :: posvar |
---|
| 101 | INTEGER, DIMENSION(:) :: point |
---|
| 102 | C |
---|
| 103 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 104 | INTEGER :: dimensio ! DIMENSION of the variable |
---|
| 105 | INTEGER :: i |
---|
| 106 | C |
---|
| 107 | C |
---|
| 108 | C Begin |
---|
| 109 | C |
---|
| 110 | dimensio = Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim |
---|
| 111 | C |
---|
| 112 | if (.not.associated(Agrif_Mygrid % tabvars(tabvarsindic) |
---|
| 113 | & %var % posvar)) then |
---|
| 114 | Allocate( |
---|
| 115 | & Agrif_Mygrid % tabvars(tabvarsindic)%var % posvar(dimensio)) |
---|
| 116 | endif |
---|
| 117 | |
---|
| 118 | do i = 1 , dimensio |
---|
| 119 | Agrif_Mygrid % tabvars(tabvarsindic) %var % posvar(i) |
---|
| 120 | & = posvar(i) |
---|
| 121 | Agrif_Mygrid % tabvars(tabvarsindic) %var % point(i) |
---|
| 122 | & = point(i) |
---|
| 123 | enddo |
---|
| 124 | C |
---|
| 125 | C |
---|
| 126 | End Subroutine Agrif_Set_type |
---|
| 127 | C |
---|
| 128 | C |
---|
| 129 | C ************************************************************************** |
---|
| 130 | CCC Subroutine Agrif_Set_parent_int |
---|
| 131 | C ************************************************************************** |
---|
| 132 | C |
---|
| 133 | Subroutine Agrif_Set_parent_int(tabvarsindic,value) |
---|
| 134 | C |
---|
| 135 | CCC Description: |
---|
| 136 | CCC To set the TYPE of the variable. |
---|
| 137 | C |
---|
| 138 | C Modules used: |
---|
| 139 | C |
---|
| 140 | |
---|
| 141 | C |
---|
| 142 | C Declarations: |
---|
| 143 | C |
---|
| 144 | C |
---|
| 145 | C |
---|
| 146 | C Arguments |
---|
| 147 | C |
---|
| 148 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 149 | INTEGER :: Value |
---|
| 150 | C |
---|
| 151 | C Begin |
---|
| 152 | C |
---|
| 153 | Agrif_Curgrid % parent % tabvars(tabvarsindic) % |
---|
| 154 | & var % iarray0 = value |
---|
| 155 | C |
---|
| 156 | C |
---|
| 157 | End Subroutine Agrif_Set_parent_int |
---|
| 158 | C |
---|
| 159 | C |
---|
| 160 | C ************************************************************************** |
---|
| 161 | CCC Subroutine Agrif_Set_parent_real |
---|
| 162 | C ************************************************************************** |
---|
| 163 | C |
---|
| 164 | Subroutine Agrif_Set_parent_real(tabvarsindic,value) |
---|
| 165 | C |
---|
| 166 | CCC Description: |
---|
| 167 | CCC To set the TYPE of the variable. |
---|
| 168 | C |
---|
| 169 | C Modules used: |
---|
| 170 | C |
---|
| 171 | |
---|
| 172 | C |
---|
| 173 | C Declarations: |
---|
| 174 | C |
---|
| 175 | C |
---|
| 176 | C |
---|
| 177 | C Arguments |
---|
| 178 | C |
---|
| 179 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 180 | REAL :: Value |
---|
| 181 | C |
---|
| 182 | C Begin |
---|
| 183 | C |
---|
| 184 | Agrif_Curgrid % parent % tabvars(tabvarsindic) % |
---|
| 185 | & var % array0 = value |
---|
| 186 | C |
---|
| 187 | C |
---|
| 188 | End Subroutine Agrif_Set_parent_real |
---|
| 189 | C |
---|
| 190 | C |
---|
| 191 | C |
---|
| 192 | C ************************************************************************** |
---|
| 193 | CCC Subroutine Agrif_Set_raf |
---|
| 194 | C ************************************************************************** |
---|
| 195 | C |
---|
| 196 | Subroutine Agrif_Set_raf(tabvarsindic,tabraf) |
---|
| 197 | C |
---|
| 198 | CCC Description: |
---|
| 199 | CCC Attention tabraf est de taille trois si on ne raffine pas suivant z la |
---|
| 200 | CCC troisieme entree du tableau tabraf est 'N' |
---|
| 201 | C |
---|
| 202 | C Modules used: |
---|
| 203 | C |
---|
| 204 | |
---|
| 205 | C |
---|
| 206 | C Declarations: |
---|
| 207 | C |
---|
| 208 | C Arguments |
---|
| 209 | C |
---|
| 210 | CHARACTER(*) ,DIMENSION(:) :: tabraf |
---|
| 211 | C |
---|
| 212 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 213 | INTEGER :: dimensio ! DIMENSION of the variable |
---|
| 214 | INTEGER :: i |
---|
| 215 | C |
---|
| 216 | C |
---|
| 217 | C Begin |
---|
| 218 | C |
---|
| 219 | dimensio = Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim |
---|
| 220 | C |
---|
| 221 | if (.not.associated(Agrif_Mygrid % tabvars(tabvarsindic) |
---|
| 222 | & %var % interptab)) then |
---|
| 223 | Allocate( |
---|
| 224 | & Agrif_Mygrid % tabvars(tabvarsindic)%var% interptab(dimensio)) |
---|
| 225 | endif |
---|
| 226 | |
---|
| 227 | do i = 1 , dimensio |
---|
| 228 | Agrif_Mygrid % tabvars(tabvarsindic) %var % interptab(i) |
---|
| 229 | & = TRIM(tabraf(i)) |
---|
| 230 | enddo |
---|
| 231 | C |
---|
| 232 | End Subroutine Agrif_Set_raf |
---|
| 233 | C |
---|
| 234 | C |
---|
| 235 | C |
---|
| 236 | C ************************************************************************** |
---|
| 237 | CCC Subroutine Agrif_Set_bc |
---|
| 238 | C ************************************************************************** |
---|
| 239 | C |
---|
| 240 | Subroutine Agrif_Set_bc(tabvarsindic,point, |
---|
| 241 | & Interpolationshouldbemade) |
---|
| 242 | C |
---|
| 243 | CCC Description: |
---|
| 244 | CCC |
---|
| 245 | C |
---|
| 246 | C Modules used: |
---|
| 247 | C |
---|
| 248 | |
---|
| 249 | C |
---|
| 250 | C Declarations: |
---|
| 251 | C |
---|
| 252 | C Arguments |
---|
| 253 | C |
---|
| 254 | INTEGER, DIMENSION(2) :: point |
---|
| 255 | LOGICAL, OPTIONAL :: Interpolationshouldbemade |
---|
| 256 | C |
---|
| 257 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 258 | TYPE(Agrif_PVariable),Pointer ::tabvars |
---|
| 259 | |
---|
| 260 | |
---|
| 261 | C |
---|
| 262 | C |
---|
| 263 | C Begin |
---|
| 264 | C |
---|
| 265 | C |
---|
| 266 | |
---|
| 267 | if (tabvarsindic <=0) then |
---|
| 268 | tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) |
---|
| 269 | else |
---|
| 270 | tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) |
---|
| 271 | endif |
---|
| 272 | |
---|
| 273 | if (Agrif_Curgrid % fixedrank .NE. 0) then |
---|
| 274 | IF (.Not.Associated(tabvars%var% interpIndex)) THEN |
---|
| 275 | Allocate(tabvars%var % interpIndex) |
---|
| 276 | tabvars%var % interpIndex = -1 |
---|
| 277 | |
---|
| 278 | Allocate(tabvars%var % oldvalues2D(2,1)) |
---|
| 279 | tabvars%var % oldvalues2D = 0. |
---|
| 280 | ENDIF |
---|
| 281 | if ( PRESENT(Interpolationshouldbemade) ) then |
---|
| 282 | tabvars%var % |
---|
| 283 | & Interpolationshouldbemade = Interpolationshouldbemade |
---|
| 284 | endif |
---|
| 285 | |
---|
| 286 | endif |
---|
| 287 | C |
---|
| 288 | tabvars%var % bcinf = point(1) |
---|
| 289 | tabvars%var % bcsup = point(2) |
---|
| 290 | C |
---|
| 291 | End Subroutine Agrif_Set_bc |
---|
| 292 | C |
---|
| 293 | C |
---|
| 294 | C ************************************************************************** |
---|
| 295 | CCC Subroutine Agrif_Set_interp |
---|
| 296 | C ************************************************************************** |
---|
| 297 | C |
---|
| 298 | Subroutine Agrif_Set_interp(tabvarsindic,interp,interp1,interp2, |
---|
| 299 | & interp3) |
---|
| 300 | C |
---|
| 301 | CCC Description: |
---|
| 302 | C |
---|
| 303 | C Declarations: |
---|
| 304 | C |
---|
| 305 | C Arguments |
---|
| 306 | C |
---|
| 307 | INTEGER, OPTIONAL :: interp,interp1,interp2,interp3 |
---|
| 308 | C |
---|
| 309 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 310 | C |
---|
| 311 | C Begin |
---|
| 312 | C |
---|
| 313 | Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp = |
---|
| 314 | & Agrif_Constant |
---|
| 315 | IF (present(interp)) THEN |
---|
| 316 | Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp = |
---|
| 317 | & interp |
---|
| 318 | ENDIF |
---|
| 319 | IF (present(interp1)) THEN |
---|
| 320 | Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp(1) = |
---|
| 321 | & interp1 |
---|
| 322 | ENDIF |
---|
| 323 | IF (present(interp2)) THEN |
---|
| 324 | Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp(2) = |
---|
| 325 | & interp2 |
---|
| 326 | ENDIF |
---|
| 327 | IF (present(interp3)) THEN |
---|
| 328 | Agrif_Mygrid % tabvars(tabvarsindic)% var % Typeinterp(3) = |
---|
| 329 | & interp3 |
---|
| 330 | ENDIF |
---|
| 331 | C |
---|
| 332 | End Subroutine Agrif_Set_interp |
---|
| 333 | C |
---|
| 334 | C ************************************************************************** |
---|
| 335 | CCC Subroutine Agrif_Set_bcinterp |
---|
| 336 | C ************************************************************************** |
---|
| 337 | C |
---|
| 338 | Subroutine Agrif_Set_bcinterp(tabvarsindic,interp,interp1, |
---|
| 339 | & interp2,interp3,interp11,interp12,interp21,interp22) |
---|
| 340 | C |
---|
| 341 | CCC Description: |
---|
| 342 | |
---|
| 343 | C |
---|
| 344 | C Modules used: |
---|
| 345 | C |
---|
| 346 | |
---|
| 347 | C |
---|
| 348 | C Declarations: |
---|
| 349 | C |
---|
| 350 | C Arguments |
---|
| 351 | C |
---|
| 352 | INTEGER, OPTIONAL :: interp,interp1,interp2,interp3 |
---|
| 353 | INTEGER, OPTIONAL :: interp11,interp12,interp21,interp22 |
---|
| 354 | C |
---|
| 355 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 356 | TYPE(Agrif_PVariable),Pointer ::tabvars |
---|
| 357 | |
---|
| 358 | |
---|
| 359 | C |
---|
| 360 | C |
---|
| 361 | C Begin |
---|
| 362 | C |
---|
| 363 | C |
---|
| 364 | |
---|
| 365 | if (tabvarsindic <=0) then |
---|
| 366 | tabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) |
---|
| 367 | else |
---|
| 368 | tabvars=>Agrif_Mygrid % tabvars(tabvarsindic) |
---|
| 369 | endif |
---|
| 370 | C |
---|
| 371 | tabvars% var % bctypeinterp = |
---|
| 372 | & Agrif_Constant |
---|
| 373 | IF (present(interp)) THEN |
---|
| 374 | tabvars% var % bctypeinterp = |
---|
| 375 | & interp |
---|
| 376 | ENDIF |
---|
| 377 | IF (present(interp1)) THEN |
---|
| 378 | tabvars% var % bctypeinterp(1:2,1) = |
---|
| 379 | & interp1 |
---|
| 380 | ENDIF |
---|
| 381 | IF (present(interp11)) THEN |
---|
| 382 | tabvars% var % bctypeinterp(1,1) = |
---|
| 383 | & interp11 |
---|
| 384 | ENDIF |
---|
| 385 | IF (present(interp12)) THEN |
---|
| 386 | tabvars% var % bctypeinterp(1,2) = |
---|
| 387 | & interp12 |
---|
| 388 | ENDIF |
---|
| 389 | IF (present(interp2)) THEN |
---|
| 390 | tabvars% var % bctypeinterp(1:2,2) = |
---|
| 391 | & interp2 |
---|
| 392 | ENDIF |
---|
| 393 | IF (present(interp21)) THEN |
---|
| 394 | tabvars% var % bctypeinterp(2,1) = |
---|
| 395 | & interp21 |
---|
| 396 | ENDIF |
---|
| 397 | IF (present(interp22)) THEN |
---|
| 398 | tabvars% var % bctypeinterp(2,2) = |
---|
| 399 | & interp22 |
---|
| 400 | ENDIF |
---|
| 401 | IF (present(interp3)) THEN |
---|
| 402 | tabvars% var % bctypeinterp(1:2,3) = |
---|
| 403 | & interp3 |
---|
| 404 | ENDIF |
---|
| 405 | C |
---|
| 406 | End Subroutine Agrif_Set_bcinterp |
---|
| 407 | C |
---|
| 408 | C |
---|
| 409 | C ************************************************************************** |
---|
| 410 | CCC Subroutine Agrif_Set_Update |
---|
| 411 | C ************************************************************************** |
---|
| 412 | C |
---|
| 413 | Subroutine Agrif_Set_Update(tabvarsindic,point) |
---|
| 414 | C |
---|
| 415 | CCC Description: |
---|
| 416 | CCC |
---|
| 417 | C |
---|
| 418 | C Modules used: |
---|
| 419 | C |
---|
| 420 | |
---|
| 421 | C |
---|
| 422 | C Declarations: |
---|
| 423 | C |
---|
| 424 | C Arguments |
---|
| 425 | C |
---|
| 426 | INTEGER, DIMENSION(2) :: point |
---|
| 427 | C |
---|
| 428 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 429 | C |
---|
| 430 | C |
---|
| 431 | C Begin |
---|
| 432 | C |
---|
| 433 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf = point(1) |
---|
| 434 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup = point(2) |
---|
| 435 | C |
---|
| 436 | End Subroutine Agrif_Set_Update |
---|
| 437 | C |
---|
| 438 | C |
---|
| 439 | C |
---|
| 440 | C ************************************************************************** |
---|
| 441 | CCC Subroutine Agrif_Set_UpdateType |
---|
| 442 | C ************************************************************************** |
---|
| 443 | C |
---|
| 444 | Subroutine Agrif_Set_UpdateType(tabvarsindic, |
---|
| 445 | & update,update1,update2, |
---|
| 446 | & update3,update4,update5) |
---|
| 447 | C |
---|
| 448 | CCC Description: |
---|
| 449 | |
---|
| 450 | C |
---|
| 451 | C Modules used: |
---|
| 452 | C |
---|
| 453 | |
---|
| 454 | C |
---|
| 455 | C Declarations: |
---|
| 456 | C |
---|
| 457 | C Arguments |
---|
| 458 | C |
---|
| 459 | INTEGER, OPTIONAL :: update, update1, |
---|
| 460 | & update2, update3,update4,update5 |
---|
| 461 | C |
---|
| 462 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 463 | C |
---|
| 464 | C |
---|
| 465 | C Begin |
---|
| 466 | C |
---|
| 467 | Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate = |
---|
| 468 | & Agrif_Update_Copy |
---|
| 469 | |
---|
| 470 | IF (present(update)) THEN |
---|
| 471 | Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate = |
---|
| 472 | & update |
---|
| 473 | ENDIF |
---|
| 474 | IF (present(update1)) THEN |
---|
| 475 | Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(1) = |
---|
| 476 | & update1 |
---|
| 477 | ENDIF |
---|
| 478 | IF (present(update2)) THEN |
---|
| 479 | Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(2) = |
---|
| 480 | & update2 |
---|
| 481 | ENDIF |
---|
| 482 | IF (present(update3)) THEN |
---|
| 483 | Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(3) = |
---|
| 484 | & update3 |
---|
| 485 | ENDIF |
---|
| 486 | IF (present(update4)) THEN |
---|
| 487 | Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(4) = |
---|
| 488 | & update4 |
---|
| 489 | ENDIF |
---|
| 490 | IF (present(update5)) THEN |
---|
| 491 | Agrif_Mygrid % tabvars(tabvarsindic)% var % typeupdate(5) = |
---|
| 492 | & update5 |
---|
| 493 | ENDIF |
---|
| 494 | C |
---|
| 495 | End Subroutine Agrif_Set_UpdateType |
---|
| 496 | C |
---|
| 497 | C |
---|
| 498 | C ************************************************************************** |
---|
| 499 | CCC Subroutine Agrif_Set_restore |
---|
| 500 | C ************************************************************************** |
---|
| 501 | C |
---|
| 502 | Subroutine Agrif_Set_restore(tabvarsindic) |
---|
| 503 | C |
---|
| 504 | CCC Description: |
---|
| 505 | CCC |
---|
| 506 | C |
---|
| 507 | C Modules used: |
---|
| 508 | C |
---|
| 509 | |
---|
| 510 | C |
---|
| 511 | C Declarations: |
---|
| 512 | C |
---|
| 513 | C Arguments |
---|
| 514 | C |
---|
| 515 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 516 | C |
---|
| 517 | C Begin |
---|
| 518 | C |
---|
| 519 | C |
---|
| 520 | Agrif_Mygrid%tabvars(tabvarsindic)%var % restaure = .TRUE. |
---|
| 521 | C |
---|
| 522 | End Subroutine Agrif_Set_restore |
---|
| 523 | C |
---|
| 524 | C |
---|
| 525 | C ************************************************************************** |
---|
| 526 | CCC Subroutine Agrif_Init_variable0d |
---|
| 527 | C ************************************************************************** |
---|
| 528 | Subroutine Agrif_Init_variable0d(tabvarsindic0,tabvarsindic, |
---|
| 529 | & procname) |
---|
| 530 | |
---|
| 531 | INTEGER :: tabvarsindic0 ! indice of the variable in tabvars |
---|
| 532 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 533 | External :: procname |
---|
| 534 | Optional :: procname |
---|
| 535 | C |
---|
| 536 | if (Agrif_Root()) Return |
---|
| 537 | C |
---|
| 538 | if (present(procname)) then |
---|
| 539 | CALL Agrif_Interp_variable(tabvarsindic0,tabvarsindic,procname) |
---|
| 540 | CALL Agrif_Bc_variable(tabvarsindic0,tabvarsindic,1.,procname) |
---|
| 541 | else |
---|
| 542 | CALL Agrif_Interp_variable(tabvarsindic0,tabvarsindic) |
---|
| 543 | CALL Agrif_Bc_variable(tabvarsindic0,tabvarsindic,1.) |
---|
| 544 | endif |
---|
| 545 | |
---|
| 546 | End Subroutine Agrif_Init_variable0d |
---|
| 547 | C |
---|
| 548 | C |
---|
| 549 | C ************************************************************************** |
---|
| 550 | CCC Subroutine Agrif_Init_variable1d |
---|
| 551 | C ************************************************************************** |
---|
| 552 | Subroutine Agrif_Init_variable1d(q,tabvarsindic,procname) |
---|
| 553 | |
---|
| 554 | REAL, DIMENSION(:) :: q |
---|
| 555 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 556 | External :: procname |
---|
| 557 | Optional :: procname |
---|
| 558 | |
---|
| 559 | C |
---|
| 560 | if (Agrif_Root()) Return |
---|
| 561 | C |
---|
| 562 | if (present(procname)) then |
---|
| 563 | CALL Agrif_Interp_variable(q,tabvarsindic,procname) |
---|
| 564 | CALL Agrif_Bc_variable(q,tabvarsindic,1.,procname) |
---|
| 565 | else |
---|
| 566 | CALL Agrif_Interp_variable(q,tabvarsindic) |
---|
| 567 | CALL Agrif_Bc_variable(q,tabvarsindic,1.) |
---|
| 568 | endif |
---|
| 569 | |
---|
| 570 | End Subroutine Agrif_Init_variable1d |
---|
| 571 | C |
---|
| 572 | C ************************************************************************** |
---|
| 573 | CCC Subroutine Agrif_Init_variable2d |
---|
| 574 | C ************************************************************************** |
---|
| 575 | Subroutine Agrif_Init_variable2d(q,tabvarsindic,procname) |
---|
| 576 | |
---|
| 577 | REAL, DIMENSION(:,:) :: q |
---|
| 578 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 579 | External :: procname |
---|
| 580 | Optional :: procname |
---|
| 581 | |
---|
| 582 | C |
---|
| 583 | if (Agrif_Root()) Return |
---|
| 584 | C |
---|
| 585 | if (present(procname)) then |
---|
| 586 | CALL Agrif_Interp_variable(q,tabvarsindic,procname) |
---|
| 587 | CALL Agrif_Bc_variable(q,tabvarsindic,1.,procname) |
---|
| 588 | else |
---|
| 589 | CALL Agrif_Interp_variable(q,tabvarsindic) |
---|
| 590 | CALL Agrif_Bc_variable(q,tabvarsindic,1.) |
---|
| 591 | endif |
---|
| 592 | |
---|
| 593 | |
---|
| 594 | End Subroutine Agrif_Init_variable2d |
---|
| 595 | C |
---|
| 596 | C |
---|
| 597 | C ************************************************************************** |
---|
| 598 | CCC Subroutine Agrif_Init_variable3d |
---|
| 599 | C ************************************************************************** |
---|
| 600 | Subroutine Agrif_Init_variable3d(q,tabvarsindic,procname) |
---|
| 601 | |
---|
| 602 | REAL, DIMENSION(:,:,:) :: q |
---|
| 603 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 604 | External :: procname |
---|
| 605 | Optional :: procname |
---|
| 606 | C |
---|
| 607 | if (Agrif_Root()) Return |
---|
| 608 | C |
---|
| 609 | if (present(procname)) then |
---|
| 610 | CALL Agrif_Interp_variable(q,tabvarsindic,procname) |
---|
| 611 | CALL Agrif_Bc_variable(q,tabvarsindic,1.,procname) |
---|
| 612 | else |
---|
| 613 | CALL Agrif_Interp_variable(q,tabvarsindic) |
---|
| 614 | CALL Agrif_Bc_variable(q,tabvarsindic,1.) |
---|
| 615 | endif |
---|
| 616 | |
---|
| 617 | C |
---|
| 618 | End Subroutine Agrif_Init_variable3d |
---|
| 619 | C |
---|
| 620 | C |
---|
| 621 | C ************************************************************************** |
---|
| 622 | CCC Subroutine Agrif_Init_variable4d |
---|
| 623 | C ************************************************************************** |
---|
| 624 | Subroutine Agrif_Init_variable4d(q,tabvarsindic,procname) |
---|
| 625 | |
---|
| 626 | REAL, DIMENSION(:,:,:,:) :: q |
---|
| 627 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 628 | External :: procname |
---|
| 629 | Optional :: procname |
---|
| 630 | C |
---|
| 631 | if (Agrif_Root()) Return |
---|
| 632 | C |
---|
| 633 | if (present(procname)) then |
---|
| 634 | CALL Agrif_Interp_variable(q,tabvarsindic,procname) |
---|
| 635 | CALL Agrif_Bc_variable(q,tabvarsindic,1.,procname) |
---|
| 636 | else |
---|
| 637 | CALL Agrif_Interp_variable(q,tabvarsindic) |
---|
| 638 | CALL Agrif_Bc_variable(q,tabvarsindic,1.) |
---|
| 639 | endif |
---|
| 640 | |
---|
| 641 | C |
---|
| 642 | End Subroutine Agrif_Init_variable4d |
---|
| 643 | C |
---|
| 644 | C |
---|
| 645 | C ************************************************************************** |
---|
| 646 | CCC Subroutine Agrif_Bc_variable0d |
---|
| 647 | C ************************************************************************** |
---|
| 648 | Subroutine Agrif_Bc_variable0d(tabvarsindic0,tabvarsindic, |
---|
| 649 | & calledweight,procname) |
---|
| 650 | |
---|
| 651 | INTEGER :: tabvarsindic0 ! indice of the variable in tabvars |
---|
| 652 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 653 | C |
---|
| 654 | External :: procname |
---|
| 655 | Optional :: procname |
---|
| 656 | REAL, OPTIONAL :: calledweight |
---|
| 657 | REAL :: weight |
---|
| 658 | LOGICAL :: pweight |
---|
| 659 | C |
---|
| 660 | INTEGER :: dimensio |
---|
| 661 | |
---|
| 662 | if (Agrif_Root()) Return |
---|
| 663 | C |
---|
| 664 | dimensio = Agrif_Mygrid % tabvars(tabvarsindic) % var % nbdim |
---|
| 665 | C |
---|
| 666 | if ( PRESENT(calledweight) ) then |
---|
| 667 | weight=calledweight |
---|
| 668 | pweight = .TRUE. |
---|
| 669 | else |
---|
| 670 | weight = 0. |
---|
| 671 | pweight = .FALSE. |
---|
| 672 | endif |
---|
| 673 | C |
---|
| 674 | C |
---|
| 675 | |
---|
| 676 | |
---|
| 677 | if ( dimensio .EQ. 1 ) Call Agrif_Interp_Bc_1D( |
---|
| 678 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 679 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 680 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 681 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array1, |
---|
| 682 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 683 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 684 | & weight, |
---|
| 685 | & pweight) |
---|
| 686 | C |
---|
| 687 | if ( dimensio .EQ. 2 ) then |
---|
| 688 | IF (present(procname)) THEN |
---|
| 689 | Call Agrif_Interp_Bc_2D( |
---|
| 690 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 691 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 692 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 693 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array2, |
---|
| 694 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 695 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 696 | & weight,pweight,procname) |
---|
| 697 | ELSE |
---|
| 698 | |
---|
| 699 | Call Agrif_Interp_Bc_2D( |
---|
| 700 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 701 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 702 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 703 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array2, |
---|
| 704 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 705 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 706 | & weight,pweight) |
---|
| 707 | ENDIF |
---|
| 708 | endif |
---|
| 709 | C |
---|
| 710 | if ( dimensio .EQ. 3 ) then |
---|
| 711 | IF (present(procname)) THEN |
---|
| 712 | |
---|
| 713 | Call Agrif_Interp_Bc_3D( |
---|
| 714 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 715 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 716 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 717 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array3, |
---|
| 718 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 719 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 720 | & weight,pweight,procname) |
---|
| 721 | ELSE |
---|
| 722 | Call Agrif_Interp_Bc_3D( |
---|
| 723 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 724 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 725 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 726 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array3, |
---|
| 727 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 728 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 729 | & weight,pweight) |
---|
| 730 | ENDIF |
---|
| 731 | endif |
---|
| 732 | C |
---|
| 733 | if ( dimensio .EQ. 4 ) then |
---|
| 734 | IF (present(procname)) THEN |
---|
| 735 | Call Agrif_Interp_Bc_4D( |
---|
| 736 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 737 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 738 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 739 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array4, |
---|
| 740 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 741 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 742 | & weight,pweight,procname) |
---|
| 743 | ELSE |
---|
| 744 | Call Agrif_Interp_Bc_4D( |
---|
| 745 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 746 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 747 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 748 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array4, |
---|
| 749 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 750 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 751 | & weight,pweight) |
---|
| 752 | ENDIF |
---|
| 753 | endif |
---|
| 754 | C |
---|
| 755 | if ( dimensio .EQ. 5 ) then |
---|
| 756 | IF (present(procname)) THEN |
---|
| 757 | Call Agrif_Interp_Bc_5D( |
---|
| 758 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 759 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 760 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 761 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array5, |
---|
| 762 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 763 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 764 | & weight,pweight,procname) |
---|
| 765 | ELSE |
---|
| 766 | Call Agrif_Interp_Bc_5D( |
---|
| 767 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 768 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 769 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 770 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array5, |
---|
| 771 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 772 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 773 | & weight,pweight) |
---|
| 774 | ENDIF |
---|
| 775 | endif |
---|
| 776 | C |
---|
| 777 | if ( dimensio .EQ. 6 ) Call Agrif_Interp_Bc_6D( |
---|
| 778 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % bctypeinterp, |
---|
| 779 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 780 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 781 | & Agrif_Curgrid % tabvars(tabvarsindic0) %var % array6, |
---|
| 782 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcinf, |
---|
| 783 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % bcsup, |
---|
| 784 | & weight, |
---|
| 785 | & pweight) |
---|
| 786 | C |
---|
| 787 | C |
---|
| 788 | End Subroutine Agrif_Bc_variable0d |
---|
| 789 | C |
---|
| 790 | C |
---|
| 791 | C ************************************************************************** |
---|
| 792 | CCC Subroutine Agrif_Bc_variable1d |
---|
| 793 | C ************************************************************************** |
---|
| 794 | Subroutine Agrif_Bc_variable1d(q,tabvarsindic,calledweight, |
---|
| 795 | & procname) |
---|
| 796 | |
---|
| 797 | REAL , Dimension(:) :: q |
---|
| 798 | External :: procname |
---|
| 799 | Optional :: procname |
---|
| 800 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 801 | C |
---|
| 802 | REAL, OPTIONAL :: calledweight |
---|
| 803 | REAL :: weight |
---|
| 804 | LOGICAL :: pweight |
---|
| 805 | TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars |
---|
| 806 | C |
---|
| 807 | C |
---|
| 808 | C |
---|
| 809 | If (Agrif_Root()) Return |
---|
| 810 | |
---|
| 811 | if ( PRESENT(calledweight) ) then |
---|
| 812 | weight=calledweight |
---|
| 813 | pweight = .TRUE. |
---|
| 814 | else |
---|
| 815 | weight = 0. |
---|
| 816 | pweight = .FALSE. |
---|
| 817 | endif |
---|
| 818 | |
---|
| 819 | if (tabvarsindic <=0) then |
---|
| 820 | tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) |
---|
| 821 | parenttabvars => tabvars%parent_var |
---|
| 822 | roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) |
---|
| 823 | else |
---|
| 824 | tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) |
---|
| 825 | parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) |
---|
| 826 | roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) |
---|
| 827 | endif |
---|
| 828 | |
---|
| 829 | IF (present(procname)) THEN |
---|
| 830 | Call Agrif_Interp_Bc_1D( |
---|
| 831 | & roottabvars % var % bctypeinterp, |
---|
| 832 | & parenttabvars, |
---|
| 833 | & tabvars,q, |
---|
| 834 | & tabvars % var % bcinf, |
---|
| 835 | & tabvars % var % bcsup, |
---|
| 836 | & weight,pweight,procname) |
---|
| 837 | ELSE |
---|
| 838 | Call Agrif_Interp_Bc_1D( |
---|
| 839 | & roottabvars % var % bctypeinterp, |
---|
| 840 | & parenttabvars, |
---|
| 841 | & tabvars,q, |
---|
| 842 | & tabvars % var % bcinf, |
---|
| 843 | & tabvars % var % bcsup, |
---|
| 844 | & weight,pweight) |
---|
| 845 | ENDIF |
---|
| 846 | End Subroutine Agrif_Bc_variable1d |
---|
| 847 | |
---|
| 848 | C |
---|
| 849 | C ************************************************************************** |
---|
| 850 | CCC Subroutine Agrif_Bc_variable2d |
---|
| 851 | C ************************************************************************** |
---|
| 852 | Subroutine Agrif_Bc_variable2d(q,tabvarsindic,calledweight, |
---|
| 853 | & procname) |
---|
| 854 | |
---|
| 855 | REAL , Dimension(:,:) :: q |
---|
| 856 | External :: procname |
---|
| 857 | Optional :: procname |
---|
| 858 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 859 | C |
---|
| 860 | REAL, OPTIONAL :: calledweight |
---|
| 861 | REAL :: weight |
---|
| 862 | LOGICAL :: pweight |
---|
| 863 | TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars |
---|
| 864 | C |
---|
| 865 | C |
---|
| 866 | C |
---|
| 867 | If (Agrif_Root()) Return |
---|
| 868 | |
---|
| 869 | if ( PRESENT(calledweight) ) then |
---|
| 870 | weight=calledweight |
---|
| 871 | pweight = .TRUE. |
---|
| 872 | else |
---|
| 873 | weight = 0. |
---|
| 874 | pweight = .FALSE. |
---|
| 875 | endif |
---|
| 876 | |
---|
| 877 | if (tabvarsindic <=0) then |
---|
| 878 | tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) |
---|
| 879 | parenttabvars => tabvars%parent_var |
---|
| 880 | roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) |
---|
| 881 | else |
---|
| 882 | tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) |
---|
| 883 | parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) |
---|
| 884 | roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) |
---|
| 885 | endif |
---|
| 886 | |
---|
| 887 | IF (present(procname)) THEN |
---|
| 888 | Call Agrif_Interp_Bc_2D( |
---|
| 889 | & roottabvars % var % bctypeinterp, |
---|
| 890 | & parenttabvars, |
---|
| 891 | & tabvars,q, |
---|
| 892 | & tabvars % var % bcinf, |
---|
| 893 | & tabvars % var % bcsup, |
---|
| 894 | & weight,pweight,procname) |
---|
| 895 | ELSE |
---|
| 896 | Call Agrif_Interp_Bc_2D( |
---|
| 897 | & roottabvars % var % bctypeinterp, |
---|
| 898 | & parenttabvars, |
---|
| 899 | & tabvars,q, |
---|
| 900 | & tabvars % var % bcinf, |
---|
| 901 | & tabvars % var % bcsup, |
---|
| 902 | & weight,pweight) |
---|
| 903 | ENDIF |
---|
| 904 | End Subroutine Agrif_Bc_variable2d |
---|
| 905 | |
---|
| 906 | C |
---|
| 907 | C ************************************************************************** |
---|
| 908 | CCC Subroutine Agrif_Bc_variable3d |
---|
| 909 | C ************************************************************************** |
---|
| 910 | Subroutine Agrif_Bc_variable3d(q,tabvarsindic,calledweight, |
---|
| 911 | & procname) |
---|
| 912 | |
---|
| 913 | REAL , Dimension(:,:,:) :: q |
---|
| 914 | External :: procname |
---|
| 915 | Optional :: procname |
---|
| 916 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 917 | C |
---|
| 918 | REAL, OPTIONAL :: calledweight |
---|
| 919 | REAL :: weight |
---|
| 920 | LOGICAL :: pweight |
---|
| 921 | TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars |
---|
| 922 | C |
---|
| 923 | C |
---|
| 924 | C |
---|
| 925 | If (Agrif_Root()) Return |
---|
| 926 | |
---|
| 927 | if ( PRESENT(calledweight) ) then |
---|
| 928 | weight=calledweight |
---|
| 929 | pweight = .TRUE. |
---|
| 930 | else |
---|
| 931 | weight = 0. |
---|
| 932 | pweight = .FALSE. |
---|
| 933 | endif |
---|
| 934 | |
---|
| 935 | if (tabvarsindic <=0) then |
---|
| 936 | tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) |
---|
| 937 | parenttabvars => tabvars%parent_var |
---|
| 938 | roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) |
---|
| 939 | else |
---|
| 940 | tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) |
---|
| 941 | parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) |
---|
| 942 | roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) |
---|
| 943 | endif |
---|
| 944 | |
---|
| 945 | IF (present(procname)) THEN |
---|
| 946 | Call Agrif_Interp_Bc_3D( |
---|
| 947 | & roottabvars % var % bctypeinterp, |
---|
| 948 | & parenttabvars, |
---|
| 949 | & tabvars,q, |
---|
| 950 | & tabvars % var % bcinf, |
---|
| 951 | & tabvars % var % bcsup, |
---|
| 952 | & weight,pweight,procname) |
---|
| 953 | ELSE |
---|
| 954 | Call Agrif_Interp_Bc_3D( |
---|
| 955 | & roottabvars % var % bctypeinterp, |
---|
| 956 | & parenttabvars, |
---|
| 957 | & tabvars,q, |
---|
| 958 | & tabvars % var % bcinf, |
---|
| 959 | & tabvars % var % bcsup, |
---|
| 960 | & weight,pweight) |
---|
| 961 | ENDIF |
---|
| 962 | End Subroutine Agrif_Bc_variable3d |
---|
| 963 | |
---|
| 964 | C |
---|
| 965 | C ************************************************************************** |
---|
| 966 | CCC Subroutine Agrif_Bc_variable4d |
---|
| 967 | C ************************************************************************** |
---|
| 968 | Subroutine Agrif_Bc_variable4d(q,tabvarsindic,calledweight, |
---|
| 969 | & procname) |
---|
| 970 | |
---|
| 971 | REAL , Dimension(:,:,:,:) :: q |
---|
| 972 | External :: procname |
---|
| 973 | Optional :: procname |
---|
| 974 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 975 | C |
---|
| 976 | REAL, OPTIONAL :: calledweight |
---|
| 977 | REAL :: weight |
---|
| 978 | LOGICAL :: pweight |
---|
| 979 | TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars |
---|
| 980 | C |
---|
| 981 | C |
---|
| 982 | C |
---|
| 983 | If (Agrif_Root()) Return |
---|
| 984 | |
---|
| 985 | if ( PRESENT(calledweight) ) then |
---|
| 986 | weight=calledweight |
---|
| 987 | pweight = .TRUE. |
---|
| 988 | else |
---|
| 989 | weight = 0. |
---|
| 990 | pweight = .FALSE. |
---|
| 991 | endif |
---|
| 992 | |
---|
| 993 | if (tabvarsindic <=0) then |
---|
| 994 | tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) |
---|
| 995 | parenttabvars => tabvars%parent_var |
---|
| 996 | roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) |
---|
| 997 | else |
---|
| 998 | tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) |
---|
| 999 | parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) |
---|
| 1000 | roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) |
---|
| 1001 | endif |
---|
| 1002 | |
---|
| 1003 | IF (present(procname)) THEN |
---|
| 1004 | Call Agrif_Interp_Bc_4D( |
---|
| 1005 | & roottabvars % var % bctypeinterp, |
---|
| 1006 | & parenttabvars, |
---|
| 1007 | & tabvars,q, |
---|
| 1008 | & tabvars % var % bcinf, |
---|
| 1009 | & tabvars % var % bcsup, |
---|
| 1010 | & weight,pweight,procname) |
---|
| 1011 | ELSE |
---|
| 1012 | Call Agrif_Interp_Bc_4D( |
---|
| 1013 | & roottabvars % var % bctypeinterp, |
---|
| 1014 | & parenttabvars, |
---|
| 1015 | & tabvars,q, |
---|
| 1016 | & tabvars % var % bcinf, |
---|
| 1017 | & tabvars % var % bcsup, |
---|
| 1018 | & weight,pweight) |
---|
| 1019 | ENDIF |
---|
| 1020 | End Subroutine Agrif_Bc_variable4d |
---|
| 1021 | |
---|
| 1022 | C |
---|
| 1023 | C ************************************************************************** |
---|
| 1024 | CCC Subroutine Agrif_Bc_variable5d |
---|
| 1025 | C ************************************************************************** |
---|
| 1026 | Subroutine Agrif_Bc_variable5d(q,tabvarsindic,calledweight, |
---|
| 1027 | & procname) |
---|
| 1028 | |
---|
| 1029 | REAL , Dimension(:,:,:,:,:) :: q |
---|
| 1030 | External :: procname |
---|
| 1031 | Optional :: procname |
---|
| 1032 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1033 | C |
---|
| 1034 | REAL, OPTIONAL :: calledweight |
---|
| 1035 | REAL :: weight |
---|
| 1036 | LOGICAL :: pweight |
---|
| 1037 | TYPE(Agrif_PVariable),Pointer ::tabvars,parenttabvars,roottabvars |
---|
| 1038 | C |
---|
| 1039 | C |
---|
| 1040 | C |
---|
| 1041 | If (Agrif_Root()) Return |
---|
| 1042 | |
---|
| 1043 | if ( PRESENT(calledweight) ) then |
---|
| 1044 | weight=calledweight |
---|
| 1045 | pweight = .TRUE. |
---|
| 1046 | else |
---|
| 1047 | weight = 0. |
---|
| 1048 | pweight = .FALSE. |
---|
| 1049 | endif |
---|
| 1050 | |
---|
| 1051 | if (tabvarsindic <=0) then |
---|
| 1052 | tabvars => Agrif_Search_Variable(Agrif_Curgrid,-tabvarsindic) |
---|
| 1053 | parenttabvars => tabvars%parent_var |
---|
| 1054 | roottabvars => Agrif_Search_Variable(Agrif_Mygrid,-tabvarsindic) |
---|
| 1055 | else |
---|
| 1056 | tabvars=>Agrif_Curgrid % tabvars(tabvarsindic) |
---|
| 1057 | parenttabvars => Agrif_Curgrid % parent % tabvars(tabvarsindic) |
---|
| 1058 | roottabvars => Agrif_Mygrid % tabvars(tabvarsindic) |
---|
| 1059 | endif |
---|
| 1060 | |
---|
| 1061 | IF (present(procname)) THEN |
---|
| 1062 | Call Agrif_Interp_Bc_5d( |
---|
| 1063 | & roottabvars % var % bctypeinterp, |
---|
| 1064 | & parenttabvars, |
---|
| 1065 | & tabvars,q, |
---|
| 1066 | & tabvars % var % bcinf, |
---|
| 1067 | & tabvars % var % bcsup, |
---|
| 1068 | & weight,pweight,procname) |
---|
| 1069 | ELSE |
---|
| 1070 | Call Agrif_Interp_Bc_5d( |
---|
| 1071 | & roottabvars % var % bctypeinterp, |
---|
| 1072 | & parenttabvars, |
---|
| 1073 | & tabvars,q, |
---|
| 1074 | & tabvars % var % bcinf, |
---|
| 1075 | & tabvars % var % bcsup, |
---|
| 1076 | & weight,pweight) |
---|
| 1077 | ENDIF |
---|
| 1078 | End Subroutine Agrif_Bc_variable5d |
---|
| 1079 | |
---|
| 1080 | C |
---|
| 1081 | C ************************************************************************** |
---|
| 1082 | CCC Subroutine Agrif_Interp_var0D |
---|
| 1083 | C ************************************************************************** |
---|
| 1084 | C |
---|
| 1085 | Subroutine Agrif_Interp_var0d(tabvarsindic0,tabvarsindic,procname) |
---|
| 1086 | |
---|
| 1087 | INTEGER :: tabvarsindic0 ! indice of the variable in tabvars |
---|
| 1088 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1089 | INTEGER :: dimensio ! indice of the variable in tabvars |
---|
| 1090 | External :: procname |
---|
| 1091 | Optional :: procname |
---|
| 1092 | C |
---|
| 1093 | if (Agrif_Root()) Return |
---|
| 1094 | C |
---|
| 1095 | dimensio = Agrif_Mygrid % tabvars(tabvarsindic) % var % nbdim |
---|
| 1096 | C |
---|
| 1097 | if ( dimensio .EQ. 1 ) then |
---|
| 1098 | if (present(procname)) then |
---|
| 1099 | Call Agrif_Interp_1D( |
---|
| 1100 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1101 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1102 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1103 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array1 , |
---|
| 1104 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1105 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1106 | else |
---|
| 1107 | Call Agrif_Interp_1D( |
---|
| 1108 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1109 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1110 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1111 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array1 , |
---|
| 1112 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1113 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1114 | endif |
---|
| 1115 | endif |
---|
| 1116 | C |
---|
| 1117 | if ( dimensio .EQ. 2 ) then |
---|
| 1118 | if (present(procname)) then |
---|
| 1119 | Call Agrif_Interp_2D( |
---|
| 1120 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1121 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1122 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1123 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array2 , |
---|
| 1124 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1125 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1126 | else |
---|
| 1127 | Call Agrif_Interp_2D( |
---|
| 1128 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1129 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1130 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1131 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array2 , |
---|
| 1132 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1133 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1134 | endif |
---|
| 1135 | endif |
---|
| 1136 | C |
---|
| 1137 | if ( dimensio .EQ. 3 ) then |
---|
| 1138 | if (present(procname)) then |
---|
| 1139 | Call Agrif_Interp_3D( |
---|
| 1140 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1141 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1142 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1143 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array3 , |
---|
| 1144 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1145 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1146 | else |
---|
| 1147 | Call Agrif_Interp_3D( |
---|
| 1148 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1149 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1150 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1151 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array3 , |
---|
| 1152 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1153 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1154 | endif |
---|
| 1155 | endif |
---|
| 1156 | C |
---|
| 1157 | if ( dimensio .EQ. 4 ) then |
---|
| 1158 | if (present(procname)) then |
---|
| 1159 | Call Agrif_Interp_4D( |
---|
| 1160 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1161 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1162 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1163 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array4 , |
---|
| 1164 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1165 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1166 | else |
---|
| 1167 | Call Agrif_Interp_4D( |
---|
| 1168 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1169 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1170 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1171 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array4 , |
---|
| 1172 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1173 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1174 | endif |
---|
| 1175 | endif |
---|
| 1176 | C |
---|
| 1177 | if ( dimensio .EQ. 5 ) then |
---|
| 1178 | if (present(procname)) then |
---|
| 1179 | Call Agrif_Interp_5D( |
---|
| 1180 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1181 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1182 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1183 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array5 , |
---|
| 1184 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1185 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1186 | else |
---|
| 1187 | Call Agrif_Interp_5D( |
---|
| 1188 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1189 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1190 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1191 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array5 , |
---|
| 1192 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1193 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1194 | endif |
---|
| 1195 | endif |
---|
| 1196 | C |
---|
| 1197 | if ( dimensio .EQ. 6 ) then |
---|
| 1198 | if (present(procname)) then |
---|
| 1199 | Call Agrif_Interp_6D( |
---|
| 1200 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1201 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1202 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1203 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array6 , |
---|
| 1204 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1205 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1206 | else |
---|
| 1207 | Call Agrif_Interp_6D( |
---|
| 1208 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1209 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1210 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1211 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array6 , |
---|
| 1212 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1213 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1214 | endif |
---|
| 1215 | endif |
---|
| 1216 | C |
---|
| 1217 | Return |
---|
| 1218 | End Subroutine Agrif_Interp_var0d |
---|
| 1219 | C |
---|
| 1220 | C ************************************************************************** |
---|
| 1221 | CCC Subroutine Agrif_Interp_var1d |
---|
| 1222 | C ************************************************************************** |
---|
| 1223 | C |
---|
| 1224 | Subroutine Agrif_Interp_var1d(q,tabvarsindic,procname) |
---|
| 1225 | |
---|
| 1226 | REAL, DIMENSION(:) :: q |
---|
| 1227 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1228 | External :: procname |
---|
| 1229 | Optional :: procname |
---|
| 1230 | C |
---|
| 1231 | if (Agrif_Root()) Return |
---|
| 1232 | C |
---|
| 1233 | if (present(procname)) then |
---|
| 1234 | Call Agrif_Interp_1D( |
---|
| 1235 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1236 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1237 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1238 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1239 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1240 | else |
---|
| 1241 | Call Agrif_Interp_1D( |
---|
| 1242 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1243 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1244 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1245 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1246 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1247 | endif |
---|
| 1248 | Return |
---|
| 1249 | End Subroutine Agrif_Interp_var1d |
---|
| 1250 | C |
---|
| 1251 | C ************************************************************************** |
---|
| 1252 | CCC Subroutine Agrif_Interp_var2d |
---|
| 1253 | C ************************************************************************** |
---|
| 1254 | C |
---|
| 1255 | Subroutine Agrif_Interp_var2d(q,tabvarsindic,procname) |
---|
| 1256 | |
---|
| 1257 | REAL, DIMENSION(:,:) :: q |
---|
| 1258 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1259 | External :: procname |
---|
| 1260 | Optional :: procname |
---|
| 1261 | |
---|
| 1262 | C |
---|
| 1263 | if (Agrif_Root()) Return |
---|
| 1264 | C |
---|
| 1265 | if (present(procname)) then |
---|
| 1266 | Call Agrif_Interp_2D( |
---|
| 1267 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1268 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1269 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1270 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1271 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1272 | else |
---|
| 1273 | Call Agrif_Interp_2D( |
---|
| 1274 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1275 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1276 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1277 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1278 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1279 | endif |
---|
| 1280 | Return |
---|
| 1281 | End Subroutine Agrif_Interp_var2d |
---|
| 1282 | C |
---|
| 1283 | C ************************************************************************** |
---|
| 1284 | CCC Subroutine Agrif_Interp_var3d |
---|
| 1285 | C ************************************************************************** |
---|
| 1286 | C |
---|
| 1287 | Subroutine Agrif_Interp_var3d(q,tabvarsindic,procname) |
---|
| 1288 | |
---|
| 1289 | REAL, DIMENSION(:,:,:) :: q |
---|
| 1290 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1291 | External :: procname |
---|
| 1292 | Optional :: procname |
---|
| 1293 | |
---|
| 1294 | C |
---|
| 1295 | if (Agrif_Root()) Return |
---|
| 1296 | C |
---|
| 1297 | if (present(procname)) then |
---|
| 1298 | Call Agrif_Interp_3D( |
---|
| 1299 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1300 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1301 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1302 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1303 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1304 | else |
---|
| 1305 | Call Agrif_Interp_3D( |
---|
| 1306 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1307 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1308 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1309 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1310 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1311 | endif |
---|
| 1312 | Return |
---|
| 1313 | End Subroutine Agrif_Interp_var3d |
---|
| 1314 | C |
---|
| 1315 | C ************************************************************************** |
---|
| 1316 | CCC Subroutine Agrif_Interp_var4d |
---|
| 1317 | C ************************************************************************** |
---|
| 1318 | C |
---|
| 1319 | Subroutine Agrif_Interp_var4d(q,tabvarsindic,procname) |
---|
| 1320 | |
---|
| 1321 | REAL, DIMENSION(:,:,:,:) :: q |
---|
| 1322 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1323 | External :: procname |
---|
| 1324 | Optional :: procname |
---|
| 1325 | |
---|
| 1326 | C |
---|
| 1327 | if (Agrif_Root()) Return |
---|
| 1328 | C |
---|
| 1329 | if (present(procname)) then |
---|
| 1330 | Call Agrif_Interp_4D( |
---|
| 1331 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1332 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1333 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1334 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1335 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1336 | else |
---|
| 1337 | Call Agrif_Interp_4D( |
---|
| 1338 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1339 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1340 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1341 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1342 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1343 | endif |
---|
| 1344 | Return |
---|
| 1345 | End Subroutine Agrif_Interp_var4d |
---|
| 1346 | C |
---|
| 1347 | C ************************************************************************** |
---|
| 1348 | CCC Subroutine Agrif_Interp_var5d |
---|
| 1349 | C ************************************************************************** |
---|
| 1350 | C |
---|
| 1351 | Subroutine Agrif_Interp_var5d(q,tabvarsindic,procname) |
---|
| 1352 | |
---|
| 1353 | REAL, DIMENSION(:,:,:,:,:) :: q |
---|
| 1354 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1355 | External :: procname |
---|
| 1356 | Optional :: procname |
---|
| 1357 | |
---|
| 1358 | C |
---|
| 1359 | if (Agrif_Root()) Return |
---|
| 1360 | C |
---|
| 1361 | if (present(procname)) then |
---|
| 1362 | Call Agrif_Interp_5D( |
---|
| 1363 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1364 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1365 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1366 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1367 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim,procname) |
---|
| 1368 | else |
---|
| 1369 | Call Agrif_Interp_5D( |
---|
| 1370 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % TypeInterp, |
---|
| 1371 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1372 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1373 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % restaure, |
---|
| 1374 | & Agrif_Mygrid % tabvars(tabvarsindic) %var % nbdim) |
---|
| 1375 | endif |
---|
| 1376 | Return |
---|
| 1377 | End Subroutine Agrif_Interp_var5d |
---|
| 1378 | C |
---|
| 1379 | C ************************************************************************** |
---|
| 1380 | CCC Subroutine Agrif_update_var0d |
---|
| 1381 | C ************************************************************************** |
---|
| 1382 | C |
---|
| 1383 | Subroutine Agrif_update_var0d(tabvarsindic0,tabvarsindic, |
---|
| 1384 | & locupdate,locupdate1, |
---|
| 1385 | & locupdate2,procname) |
---|
| 1386 | |
---|
| 1387 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1388 | INTEGER :: tabvarsindic0 ! indice of the variable in tabvars |
---|
| 1389 | External :: procname |
---|
| 1390 | Optional :: procname |
---|
| 1391 | INTEGER :: dimensio |
---|
| 1392 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate |
---|
| 1393 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 |
---|
| 1394 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate2 |
---|
| 1395 | C |
---|
| 1396 | dimensio = Agrif_Mygrid % tabvars(tabvarsindic) % var % nbdim |
---|
| 1397 | C |
---|
| 1398 | if (Agrif_Root()) Return |
---|
| 1399 | |
---|
| 1400 | C |
---|
| 1401 | IF (present(locupdate)) THEN |
---|
| 1402 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:dimensio) |
---|
| 1403 | & = locupdate(1) |
---|
| 1404 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:dimensio) |
---|
| 1405 | & = locupdate(2) |
---|
| 1406 | ELSE |
---|
| 1407 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:dimensio) |
---|
| 1408 | & = -99 |
---|
| 1409 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:dimensio) |
---|
| 1410 | & = -99 |
---|
| 1411 | ENDIF |
---|
| 1412 | |
---|
| 1413 | IF (present(locupdate1)) THEN |
---|
| 1414 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1) |
---|
| 1415 | & = locupdate1(1) |
---|
| 1416 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1) |
---|
| 1417 | & = locupdate1(2) |
---|
| 1418 | ENDIF |
---|
| 1419 | |
---|
| 1420 | IF (present(locupdate2)) THEN |
---|
| 1421 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2) |
---|
| 1422 | & = locupdate2(1) |
---|
| 1423 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2) |
---|
| 1424 | & = locupdate2(2) |
---|
| 1425 | ENDIF |
---|
| 1426 | |
---|
| 1427 | if ( dimensio .EQ. 1 ) then |
---|
| 1428 | IF (present(procname)) THEN |
---|
| 1429 | Call Agrif_Update_1D( |
---|
| 1430 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1431 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1432 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1433 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array1 , |
---|
| 1434 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1435 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1436 | & procname) |
---|
| 1437 | ELSE |
---|
| 1438 | Call Agrif_Update_1D( |
---|
| 1439 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1440 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1441 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1442 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array1 , |
---|
| 1443 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1444 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1445 | ENDIF |
---|
| 1446 | endif |
---|
| 1447 | if ( dimensio .EQ. 2 ) then |
---|
| 1448 | IF (present(procname)) THEN |
---|
| 1449 | Call Agrif_Update_2D( |
---|
| 1450 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1451 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1452 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1453 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array2 , |
---|
| 1454 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1455 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1456 | & procname) |
---|
| 1457 | ELSE |
---|
| 1458 | Call Agrif_Update_2D( |
---|
| 1459 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1460 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1461 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1462 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array2 , |
---|
| 1463 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1464 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1465 | ENDIF |
---|
| 1466 | endif |
---|
| 1467 | if ( dimensio .EQ. 3 ) then |
---|
| 1468 | IF (present(procname)) THEN |
---|
| 1469 | Call Agrif_Update_3D( |
---|
| 1470 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1471 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1472 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1473 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array3 , |
---|
| 1474 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1475 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1476 | & procname) |
---|
| 1477 | ELSE |
---|
| 1478 | Call Agrif_Update_3D( |
---|
| 1479 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1480 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1481 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1482 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array3 , |
---|
| 1483 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1484 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1485 | ENDIF |
---|
| 1486 | endif |
---|
| 1487 | if ( dimensio .EQ. 4 ) then |
---|
| 1488 | IF (present(procname)) THEN |
---|
| 1489 | Call Agrif_Update_4D( |
---|
| 1490 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1491 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1492 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1493 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array4 , |
---|
| 1494 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1495 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1496 | & procname) |
---|
| 1497 | ELSE |
---|
| 1498 | Call Agrif_Update_4D( |
---|
| 1499 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1500 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1501 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1502 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array4 , |
---|
| 1503 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1504 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1505 | ENDIF |
---|
| 1506 | endif |
---|
| 1507 | if ( dimensio .EQ. 5 ) then |
---|
| 1508 | IF (present(procname)) THEN |
---|
| 1509 | Call Agrif_Update_5D( |
---|
| 1510 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1511 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1512 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1513 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array5 , |
---|
| 1514 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1515 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1516 | & procname) |
---|
| 1517 | ELSE |
---|
| 1518 | Call Agrif_Update_5D( |
---|
| 1519 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1520 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1521 | & Agrif_Curgrid % tabvars(tabvarsindic), |
---|
| 1522 | & Agrif_Curgrid % tabvars(tabvarsindic0) % var % array5 , |
---|
| 1523 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1524 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1525 | ENDIF |
---|
| 1526 | endif |
---|
| 1527 | |
---|
| 1528 | Return |
---|
| 1529 | End Subroutine Agrif_update_var0d |
---|
| 1530 | C |
---|
| 1531 | C |
---|
| 1532 | C ************************************************************************** |
---|
| 1533 | CCC Subroutine Agrif_update_var1d |
---|
| 1534 | C ************************************************************************** |
---|
| 1535 | C |
---|
| 1536 | Subroutine Agrif_update_var1d(q,tabvarsindic,locupdate, |
---|
| 1537 | & locupdate1,locupdate2,procname) |
---|
| 1538 | |
---|
| 1539 | REAL, DIMENSION(:) :: q |
---|
| 1540 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1541 | External :: procname |
---|
| 1542 | Optional :: procname |
---|
| 1543 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate |
---|
| 1544 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 |
---|
| 1545 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate2 |
---|
| 1546 | C |
---|
| 1547 | if (Agrif_Root()) Return |
---|
| 1548 | C |
---|
| 1549 | IF (present(locupdate)) THEN |
---|
| 1550 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:1) |
---|
| 1551 | & = locupdate(1) |
---|
| 1552 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:1) |
---|
| 1553 | & = locupdate(2) |
---|
| 1554 | ELSE |
---|
| 1555 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:1) |
---|
| 1556 | & = -99 |
---|
| 1557 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:1) |
---|
| 1558 | & = -99 |
---|
| 1559 | ENDIF |
---|
| 1560 | |
---|
| 1561 | IF (present(locupdate1)) THEN |
---|
| 1562 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1) |
---|
| 1563 | & = locupdate1(1) |
---|
| 1564 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1) |
---|
| 1565 | & = locupdate1(2) |
---|
| 1566 | ENDIF |
---|
| 1567 | |
---|
| 1568 | IF (present(locupdate2)) THEN |
---|
| 1569 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2) |
---|
| 1570 | & = locupdate2(1) |
---|
| 1571 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2) |
---|
| 1572 | & = locupdate2(2) |
---|
| 1573 | ENDIF |
---|
| 1574 | |
---|
| 1575 | IF (present(procname)) THEN |
---|
| 1576 | Call Agrif_Update_1D( |
---|
| 1577 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1578 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1579 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1580 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1581 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1582 | & procname) |
---|
| 1583 | ELSE |
---|
| 1584 | Call Agrif_Update_1D( |
---|
| 1585 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1586 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1587 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1588 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1589 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1590 | ENDIF |
---|
| 1591 | |
---|
| 1592 | Return |
---|
| 1593 | End Subroutine Agrif_update_var1d |
---|
| 1594 | C |
---|
| 1595 | C |
---|
| 1596 | C ************************************************************************** |
---|
| 1597 | CCC Subroutine Agrif_update_var2d |
---|
| 1598 | C ************************************************************************** |
---|
| 1599 | C |
---|
| 1600 | Subroutine Agrif_update_var2d(q,tabvarsindic,locupdate, |
---|
| 1601 | & locupdate1,locupdate2,procname) |
---|
| 1602 | |
---|
| 1603 | REAL, DIMENSION(:,:) :: q |
---|
| 1604 | External :: procname |
---|
| 1605 | Optional :: procname |
---|
| 1606 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate |
---|
| 1607 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 |
---|
| 1608 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate2 |
---|
| 1609 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1610 | C |
---|
| 1611 | IF (Agrif_Root()) RETURN |
---|
| 1612 | |
---|
| 1613 | C |
---|
| 1614 | IF (present(locupdate)) THEN |
---|
| 1615 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:2) |
---|
| 1616 | & = locupdate(1) |
---|
| 1617 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:2) |
---|
| 1618 | & = locupdate(2) |
---|
| 1619 | ELSE |
---|
| 1620 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:2) |
---|
| 1621 | & = -99 |
---|
| 1622 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:2) |
---|
| 1623 | & = -99 |
---|
| 1624 | ENDIF |
---|
| 1625 | |
---|
| 1626 | IF (present(locupdate1)) THEN |
---|
| 1627 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1) |
---|
| 1628 | & = locupdate1(1) |
---|
| 1629 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1) |
---|
| 1630 | & = locupdate1(2) |
---|
| 1631 | ENDIF |
---|
| 1632 | |
---|
| 1633 | IF (present(locupdate2)) THEN |
---|
| 1634 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2) |
---|
| 1635 | & = locupdate2(1) |
---|
| 1636 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2) |
---|
| 1637 | & = locupdate2(2) |
---|
| 1638 | ENDIF |
---|
| 1639 | |
---|
| 1640 | IF (present(procname)) THEN |
---|
| 1641 | Call Agrif_Update_2D( |
---|
| 1642 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1643 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1644 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1645 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1646 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1647 | & procname) |
---|
| 1648 | ELSE |
---|
| 1649 | Call Agrif_Update_2D( |
---|
| 1650 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1651 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1652 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1653 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1654 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1655 | ENDIF |
---|
| 1656 | |
---|
| 1657 | Return |
---|
| 1658 | End Subroutine Agrif_update_var2d |
---|
| 1659 | C |
---|
| 1660 | C |
---|
| 1661 | C ************************************************************************** |
---|
| 1662 | CCC Subroutine Agrif_update_var3d |
---|
| 1663 | C ************************************************************************** |
---|
| 1664 | C |
---|
| 1665 | Subroutine Agrif_update_var3d(q,tabvarsindic,locupdate, |
---|
| 1666 | & locupdate1,locupdate2,procname) |
---|
| 1667 | |
---|
| 1668 | REAL, DIMENSION(:,:,:) :: q |
---|
| 1669 | External :: procname |
---|
| 1670 | Optional :: procname |
---|
| 1671 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate |
---|
| 1672 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 |
---|
| 1673 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate2 |
---|
| 1674 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1675 | C |
---|
| 1676 | IF (Agrif_Root()) RETURN |
---|
| 1677 | C |
---|
| 1678 | |
---|
| 1679 | IF (present(locupdate)) THEN |
---|
| 1680 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:3) |
---|
| 1681 | & = locupdate(1) |
---|
| 1682 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:3) |
---|
| 1683 | & = locupdate(2) |
---|
| 1684 | ELSE |
---|
| 1685 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:3) |
---|
| 1686 | & = -99 |
---|
| 1687 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:3) |
---|
| 1688 | & = -99 |
---|
| 1689 | ENDIF |
---|
| 1690 | |
---|
| 1691 | IF (present(locupdate1)) THEN |
---|
| 1692 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1) |
---|
| 1693 | & = locupdate1(1) |
---|
| 1694 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1) |
---|
| 1695 | & = locupdate1(2) |
---|
| 1696 | ENDIF |
---|
| 1697 | |
---|
| 1698 | IF (present(locupdate2)) THEN |
---|
| 1699 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2) |
---|
| 1700 | & = locupdate2(1) |
---|
| 1701 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2) |
---|
| 1702 | & = locupdate2(2) |
---|
| 1703 | ENDIF |
---|
| 1704 | |
---|
| 1705 | IF (present(procname)) THEN |
---|
| 1706 | Call Agrif_Update_3D( |
---|
| 1707 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1708 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1709 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1710 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1711 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1712 | & procname) |
---|
| 1713 | ELSE |
---|
| 1714 | Call Agrif_Update_3D( |
---|
| 1715 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1716 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1717 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1718 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1719 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1720 | ENDIF |
---|
| 1721 | |
---|
| 1722 | Return |
---|
| 1723 | End Subroutine Agrif_update_var3d |
---|
| 1724 | C |
---|
| 1725 | C |
---|
| 1726 | C ************************************************************************** |
---|
| 1727 | CCC Subroutine Agrif_update_var4d |
---|
| 1728 | C ************************************************************************** |
---|
| 1729 | C |
---|
| 1730 | Subroutine Agrif_update_var4d(q,tabvarsindic,locupdate, |
---|
| 1731 | & locupdate1,locupdate2,procname) |
---|
| 1732 | |
---|
| 1733 | REAL, DIMENSION(:,:,:,:) :: q |
---|
| 1734 | External :: procname |
---|
| 1735 | Optional :: procname |
---|
| 1736 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate |
---|
| 1737 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 |
---|
| 1738 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate2 |
---|
| 1739 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1740 | C |
---|
| 1741 | IF (Agrif_Root()) RETURN |
---|
| 1742 | C |
---|
| 1743 | IF (present(locupdate)) THEN |
---|
| 1744 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:4) |
---|
| 1745 | & = locupdate(1) |
---|
| 1746 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:4) |
---|
| 1747 | & = locupdate(2) |
---|
| 1748 | ELSE |
---|
| 1749 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:4) |
---|
| 1750 | & = -99 |
---|
| 1751 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:4) |
---|
| 1752 | & = -99 |
---|
| 1753 | ENDIF |
---|
| 1754 | |
---|
| 1755 | IF (present(locupdate1)) THEN |
---|
| 1756 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1) |
---|
| 1757 | & = locupdate1(1) |
---|
| 1758 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1) |
---|
| 1759 | & = locupdate1(2) |
---|
| 1760 | ENDIF |
---|
| 1761 | |
---|
| 1762 | IF (present(locupdate2)) THEN |
---|
| 1763 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2) |
---|
| 1764 | & = locupdate2(1) |
---|
| 1765 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2) |
---|
| 1766 | & = locupdate2(2) |
---|
| 1767 | ENDIF |
---|
| 1768 | |
---|
| 1769 | IF (present(procname)) THEN |
---|
| 1770 | Call Agrif_Update_4D( |
---|
| 1771 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1772 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1773 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1774 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1775 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1776 | & procname) |
---|
| 1777 | ELSE |
---|
| 1778 | Call Agrif_Update_4D( |
---|
| 1779 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1780 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1781 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1782 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1783 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1784 | ENDIF |
---|
| 1785 | |
---|
| 1786 | Return |
---|
| 1787 | End Subroutine Agrif_update_var4d |
---|
| 1788 | C |
---|
| 1789 | C |
---|
| 1790 | C ************************************************************************** |
---|
| 1791 | CCC Subroutine Agrif_update_var5d |
---|
| 1792 | C ************************************************************************** |
---|
| 1793 | C |
---|
| 1794 | Subroutine Agrif_update_var5d(q,tabvarsindic,locupdate, |
---|
| 1795 | & locupdate1,locupdate2,procname) |
---|
| 1796 | |
---|
| 1797 | REAL, DIMENSION(:,:,:,:,:) :: q |
---|
| 1798 | External :: procname |
---|
| 1799 | Optional :: procname |
---|
| 1800 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate |
---|
| 1801 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate1 |
---|
| 1802 | INTEGER, DIMENSION(2), OPTIONAL :: locupdate2 |
---|
| 1803 | INTEGER :: tabvarsindic ! indice of the variable in tabvars |
---|
| 1804 | C |
---|
| 1805 | IF (Agrif_Root()) RETURN |
---|
| 1806 | C |
---|
| 1807 | IF (present(locupdate)) THEN |
---|
| 1808 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:5) |
---|
| 1809 | & = locupdate(1) |
---|
| 1810 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:5) |
---|
| 1811 | & = locupdate(2) |
---|
| 1812 | ELSE |
---|
| 1813 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1:5) |
---|
| 1814 | & = -99 |
---|
| 1815 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1:5) |
---|
| 1816 | & = -99 |
---|
| 1817 | ENDIF |
---|
| 1818 | |
---|
| 1819 | IF (present(locupdate1)) THEN |
---|
| 1820 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(1) |
---|
| 1821 | & = locupdate1(1) |
---|
| 1822 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(1) |
---|
| 1823 | & = locupdate1(2) |
---|
| 1824 | ENDIF |
---|
| 1825 | |
---|
| 1826 | IF (present(locupdate2)) THEN |
---|
| 1827 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updateinf(2) |
---|
| 1828 | & = locupdate2(1) |
---|
| 1829 | Agrif_Curgrid%tabvars(tabvarsindic)%var % updatesup(2) |
---|
| 1830 | & = locupdate2(2) |
---|
| 1831 | ENDIF |
---|
| 1832 | |
---|
| 1833 | IF (present(procname)) THEN |
---|
| 1834 | Call Agrif_Update_5D( |
---|
| 1835 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1836 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1837 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1838 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1839 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup, |
---|
| 1840 | & procname) |
---|
| 1841 | ELSE |
---|
| 1842 | Call Agrif_Update_5D( |
---|
| 1843 | & Agrif_Mygrid % tabvars(tabvarsindic) % var % typeupdate, |
---|
| 1844 | & Agrif_Curgrid % parent % tabvars(tabvarsindic), |
---|
| 1845 | & Agrif_Curgrid % tabvars(tabvarsindic),q, |
---|
| 1846 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updateinf, |
---|
| 1847 | & Agrif_Curgrid % tabvars(tabvarsindic) % var % updatesup) |
---|
| 1848 | ENDIF |
---|
| 1849 | |
---|
| 1850 | Return |
---|
| 1851 | End Subroutine Agrif_update_var5d |
---|
| 1852 | |
---|
| 1853 | Subroutine Agrif_Declare_Flux(fluxname,profilename) |
---|
| 1854 | character*(*) :: fluxname, profilename |
---|
| 1855 | Type(Agrif_Flux), pointer :: newflux |
---|
| 1856 | Type(Agrif_Profile), pointer :: parcours |
---|
| 1857 | logical :: foundprofile |
---|
| 1858 | integer :: i,j,n |
---|
| 1859 | |
---|
| 1860 | foundprofile = .FALSE. |
---|
| 1861 | parcours => Agrif_Myprofiles |
---|
| 1862 | |
---|
| 1863 | Do While (Associated(parcours)) |
---|
| 1864 | IF (parcours % profilename == profilename) THEN |
---|
| 1865 | foundprofile = .TRUE. |
---|
| 1866 | EXIT |
---|
| 1867 | ENDIF |
---|
| 1868 | parcours => parcours%nextprofile |
---|
| 1869 | End Do |
---|
| 1870 | |
---|
| 1871 | IF (.NOT.foundprofile) THEN |
---|
| 1872 | write(*,*) 'The profile ''' |
---|
| 1873 | & //TRIM(profilename)//''' has not been declared' |
---|
| 1874 | stop |
---|
| 1875 | ENDIF |
---|
| 1876 | |
---|
| 1877 | Allocate(Newflux) |
---|
| 1878 | |
---|
| 1879 | Newflux % fluxname = fluxname |
---|
| 1880 | |
---|
| 1881 | Newflux % profile => parcours |
---|
| 1882 | |
---|
| 1883 | Newflux % nextflux => Agrif_Curgrid % fluxes |
---|
| 1884 | |
---|
| 1885 | Agrif_Curgrid % fluxes => Newflux |
---|
| 1886 | |
---|
| 1887 | End Subroutine Agrif_Declare_Flux |
---|
| 1888 | |
---|
| 1889 | Subroutine Agrif_Save_Flux(fluxname, fluxtab) |
---|
| 1890 | character*(*) :: fluxname |
---|
| 1891 | REAL, DIMENSION(:,:) :: fluxtab |
---|
| 1892 | |
---|
| 1893 | |
---|
| 1894 | Type(Agrif_Flux), pointer :: Flux |
---|
| 1895 | |
---|
| 1896 | Type(Agrif_pgrid), pointer :: parcours_child |
---|
| 1897 | |
---|
| 1898 | Type(Agrif_grid), Pointer :: currentgrid,oldcurgrid |
---|
| 1899 | |
---|
| 1900 | IF (.Not.Agrif_Root()) THEN |
---|
| 1901 | Flux => Agrif_Search_Flux(fluxname) |
---|
| 1902 | |
---|
| 1903 | IF (.NOT.Flux%fluxallocated) THEN |
---|
| 1904 | CALL Agrif_AllocateFlux(Flux,fluxtab) |
---|
| 1905 | ENDIF |
---|
| 1906 | |
---|
| 1907 | Call Agrif_Save_Fluxtab(Flux,fluxtab) |
---|
| 1908 | |
---|
| 1909 | ENDIF |
---|
| 1910 | |
---|
| 1911 | oldcurgrid=> Agrif_Curgrid |
---|
| 1912 | |
---|
| 1913 | parcours_child => Agrif_Curgrid%child_grids |
---|
| 1914 | |
---|
| 1915 | Do While (Associated(parcours_child)) |
---|
| 1916 | currentgrid => parcours_child%gr |
---|
| 1917 | Agrif_Curgrid => parcours_child%gr |
---|
| 1918 | Flux => Agrif_Search_Flux(fluxname) |
---|
| 1919 | IF (.NOT.Flux%fluxallocated) THEN |
---|
| 1920 | CALL Agrif_AllocateFlux(Flux,fluxtab) |
---|
| 1921 | ENDIF |
---|
| 1922 | Call Agrif_Save_Fluxtab_child(Flux,fluxtab) |
---|
| 1923 | parcours_child=> parcours_child%next |
---|
| 1924 | End Do |
---|
| 1925 | |
---|
| 1926 | Agrif_Curgrid=>oldcurgrid |
---|
| 1927 | |
---|
| 1928 | End Subroutine Agrif_Save_Flux |
---|
| 1929 | |
---|
| 1930 | Subroutine Agrif_Cancel_Flux(fluxname) |
---|
| 1931 | character*(*) :: fluxname |
---|
| 1932 | |
---|
| 1933 | Type(Agrif_Flux), pointer :: Flux |
---|
| 1934 | |
---|
| 1935 | Flux => Agrif_Search_Flux(fluxname) |
---|
| 1936 | |
---|
| 1937 | IF (Flux%FluxAllocated) Call Agrif_Cancel_Fluxarray(Flux) |
---|
| 1938 | |
---|
| 1939 | End Subroutine Agrif_Cancel_Flux |
---|
| 1940 | |
---|
| 1941 | Subroutine Agrif_Flux_Correction(fluxname, procname) |
---|
| 1942 | character*(*) :: fluxname |
---|
| 1943 | external :: procname |
---|
| 1944 | |
---|
| 1945 | Type(Agrif_Flux), pointer :: Flux |
---|
| 1946 | |
---|
| 1947 | Flux => Agrif_Search_Flux(fluxname) |
---|
| 1948 | |
---|
| 1949 | Call Agrif_FluxCorrect(Flux, procname) |
---|
| 1950 | |
---|
| 1951 | |
---|
| 1952 | End Subroutine Agrif_Flux_Correction |
---|
| 1953 | |
---|
| 1954 | Subroutine Agrif_Declare_Variable(posvar,firstpoint, |
---|
| 1955 | & raf,lb,ub,varid) |
---|
| 1956 | character*(80) :: variablename |
---|
| 1957 | Type(Agrif_List_Variables), Pointer :: newvariable,newvariablep |
---|
| 1958 | INTEGER, DIMENSION(:) :: posvar |
---|
| 1959 | INTEGER, DIMENSION(:) :: lb,ub |
---|
| 1960 | INTEGER, DIMENSION(:) :: firstpoint |
---|
| 1961 | CHARACTER(*) ,DIMENSION(:) :: raf |
---|
| 1962 | TYPE(Agrif_Pvariable), Pointer :: parent_var,root_var |
---|
| 1963 | INTEGER :: dimensio |
---|
| 1964 | INTEGER :: varid |
---|
| 1965 | |
---|
| 1966 | if (agrif_root()) return |
---|
| 1967 | |
---|
| 1968 | dimensio = SIZE(posvar) |
---|
| 1969 | C |
---|
| 1970 | C |
---|
| 1971 | Allocate(newvariable) |
---|
| 1972 | Allocate(newvariable%pvar) |
---|
| 1973 | Allocate(newvariable%pvar%var) |
---|
| 1974 | Allocate(newvariable%pvar%var%posvar(dimensio)) |
---|
| 1975 | Allocate(newvariable%pvar%var%interptab(dimensio)) |
---|
| 1976 | newvariable%pvar%var%variablename = variablename |
---|
| 1977 | newvariable%pvar%var%interptab = raf |
---|
| 1978 | newvariable%pvar%var%nbdim = dimensio |
---|
| 1979 | newvariable%pvar%var%posvar = posvar |
---|
| 1980 | newvariable%pvar%var%point(1:dimensio) = firstpoint |
---|
| 1981 | newvariable%pvar%var%lb(1:dimensio) = lb(1:dimensio) |
---|
| 1982 | newvariable%pvar%var%ub(1:dimensio) = ub(1:dimensio) |
---|
| 1983 | |
---|
| 1984 | newvariable % nextvariable => Agrif_Curgrid%variables |
---|
| 1985 | |
---|
| 1986 | Agrif_Curgrid%variables => newvariable |
---|
| 1987 | Agrif_Curgrid%Nbvariables = Agrif_Curgrid%Nbvariables + 1 |
---|
| 1988 | |
---|
| 1989 | varid = -Agrif_Curgrid%Nbvariables |
---|
| 1990 | |
---|
| 1991 | if (agrif_curgrid%parent%nbvariables < agrif_curgrid%nbvariables) |
---|
| 1992 | & then |
---|
| 1993 | Allocate(newvariablep) |
---|
| 1994 | Allocate(newvariablep%pvar) |
---|
| 1995 | Allocate(newvariablep%pvar%var) |
---|
| 1996 | Allocate(newvariablep%pvar%var%posvar(dimensio)) |
---|
| 1997 | Allocate(newvariablep%pvar%var%interptab(dimensio)) |
---|
| 1998 | newvariablep%pvar%var%variablename = variablename |
---|
| 1999 | newvariablep%pvar%var%interptab = raf |
---|
| 2000 | newvariablep%pvar%var%nbdim = dimensio |
---|
| 2001 | newvariablep%pvar%var%posvar = posvar |
---|
| 2002 | newvariablep%pvar%var%point(1:dimensio) = firstpoint |
---|
| 2003 | |
---|
| 2004 | newvariablep % nextvariable => Agrif_Curgrid%parent%variables |
---|
| 2005 | |
---|
| 2006 | Agrif_Curgrid%parent%variables => newvariablep |
---|
| 2007 | |
---|
| 2008 | Agrif_Curgrid%parent%Nbvariables = |
---|
| 2009 | & Agrif_Curgrid%parent%Nbvariables + 1 |
---|
| 2010 | parent_var=>newvariablep%pvar |
---|
| 2011 | else |
---|
| 2012 | parent_var=>Agrif_Search_Variable |
---|
| 2013 | & (Agrif_Curgrid%parent,Agrif_Curgrid%nbvariables) |
---|
| 2014 | endif |
---|
| 2015 | |
---|
| 2016 | newvariable%pvar%parent_var=>parent_var |
---|
| 2017 | |
---|
| 2018 | root_var=>Agrif_Search_Variable |
---|
| 2019 | & (Agrif_Mygrid,Agrif_Curgrid%nbvariables) |
---|
| 2020 | |
---|
| 2021 | newvariable%pvar%var%root_var=>root_var%var |
---|
| 2022 | |
---|
| 2023 | |
---|
| 2024 | End Subroutine Agrif_Declare_Variable |
---|
| 2025 | |
---|
| 2026 | FUNCTION Agrif_Search_Variable(grid,varid) |
---|
| 2027 | integer :: varid |
---|
| 2028 | Type(Agrif_Pvariable), Pointer :: Agrif_Search_variable |
---|
| 2029 | Type(Agrif_grid), Pointer :: grid |
---|
| 2030 | |
---|
| 2031 | Type(Agrif_List_Variables), pointer :: parcours |
---|
| 2032 | Logical :: foundvariable |
---|
| 2033 | integer nb |
---|
| 2034 | |
---|
| 2035 | foundvariable = .FALSE. |
---|
| 2036 | parcours => grid%variables |
---|
| 2037 | |
---|
| 2038 | do nb=1,varid-1 |
---|
| 2039 | parcours => parcours%nextvariable |
---|
| 2040 | End Do |
---|
| 2041 | |
---|
| 2042 | Agrif_Search_variable => parcours%pvar |
---|
| 2043 | |
---|
| 2044 | |
---|
| 2045 | End Function Agrif_Search_variable |
---|
| 2046 | |
---|
| 2047 | Subroutine Agrif_Declare_Profile_flux(profilename,posvar, |
---|
| 2048 | & firstpoint,raf) |
---|
| 2049 | character*(*) :: profilename |
---|
| 2050 | Type(Agrif_Profile), Pointer :: newprofile |
---|
| 2051 | INTEGER, DIMENSION(:) :: posvar |
---|
| 2052 | INTEGER, DIMENSION(:) :: firstpoint |
---|
| 2053 | CHARACTER(*) ,DIMENSION(:) :: raf |
---|
| 2054 | INTEGER :: dimensio |
---|
| 2055 | |
---|
| 2056 | dimensio = SIZE(posvar) |
---|
| 2057 | C |
---|
| 2058 | C |
---|
| 2059 | Allocate(newprofile) |
---|
| 2060 | Allocate(newprofile%posvar(dimensio)) |
---|
| 2061 | Allocate(newprofile%interptab(dimensio)) |
---|
| 2062 | newprofile%profilename = profilename |
---|
| 2063 | newprofile%interptab = raf |
---|
| 2064 | newprofile%nbdim = dimensio |
---|
| 2065 | newprofile%posvar = posvar |
---|
| 2066 | newprofile%point(1:dimensio) = firstpoint |
---|
| 2067 | |
---|
| 2068 | newprofile % nextprofile => Agrif_myprofiles |
---|
| 2069 | |
---|
| 2070 | Agrif_myprofiles => newprofile |
---|
| 2071 | |
---|
| 2072 | End Subroutine Agrif_Declare_Profile_flux |
---|
| 2073 | |
---|
| 2074 | C |
---|
| 2075 | End module Agrif_bcfunction |
---|