Changeset 11205
- Timestamp:
- 2019-07-02T14:25:46+02:00 (5 years ago)
- Location:
- NEMO/branches/2019/dev_r10973_AGRIF-01_jchanut_small_jpi_jpj/src
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/branches/2019/dev_r10973_AGRIF-01_jchanut_small_jpi_jpj/src/NST/agrif_oce.F90
r10425 r11205 46 46 47 47 ! Barotropic arrays used to store open boundary data during time-splitting loop: 48 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubdy_w, vbdy_w, hbdy_w 49 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubdy_e, vbdy_e, hbdy_e 50 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubdy_n, vbdy_n, hbdy_n 51 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubdy_s, vbdy_s, hbdy_s 48 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubdy, vbdy, hbdy 52 49 53 50 … … 91 88 & tabspongedone_v (jpi,jpj), STAT = ierr(1) ) 92 89 93 ALLOCATE( ubdy_w(nbghostcells,jpj), vbdy_w(nbghostcells,jpj), hbdy_w(nbghostcells,jpj), & 94 & ubdy_e(nbghostcells,jpj), vbdy_e(nbghostcells,jpj), hbdy_e(nbghostcells,jpj), & 95 & ubdy_n(jpi,nbghostcells), vbdy_n(jpi,nbghostcells), hbdy_n(jpi,nbghostcells), & 96 & ubdy_s(jpi,nbghostcells), vbdy_s(jpi,nbghostcells), hbdy_s(jpi,nbghostcells), STAT = ierr(2) ) 90 ALLOCATE( ubdy(jpi,jpj), vbdy(jpi,jpj), hbdy(jpi,jpj), STAT = ierr(2) ) 97 91 98 92 agrif_oce_alloc = MAXVAL(ierr) -
NEMO/branches/2019/dev_r10973_AGRIF-01_jchanut_small_jpi_jpj/src/NST/agrif_oce_interp.F90
r10068 r11205 401 401 !! 402 402 INTEGER :: ji, jj 403 INTEGER :: istart, iend, jstart, jend 403 404 !!---------------------------------------------------------------------- 404 405 ! 405 406 IF( Agrif_Root() ) RETURN 406 407 ! 407 IF((nbondi == -1).OR.(nbondi == 2)) THEN 408 !--- West ---! 409 istart = 2 410 iend = nbghostcells+1 411 DO ji = mi0(istart), mi1(iend) 408 412 DO jj=1,jpj 409 va_e( 2:nbghostcells+1,jj) = vbdy_w(1:nbghostcells,jj) * hvr_e(2:nbghostcells+1,jj)413 va_e(ji,jj) = vbdy(ji,jj) * hvr_e(ji,jj) 410 414 ! Specified fluxes: 411 ua_e(2:nbghostcells+1,jj) = ubdy_w(1:nbghostcells,jj) * hur_e(2:nbghostcells+1,jj) 412 ! Characteristics method (only if ghostcells=1): 413 !alt ua_e(2,jj) = 0.5_wp * ( ubdy_w(jj) * hur_e(2,jj) + ua_e(3,jj) & 414 !alt & - sqrt(grav * hur_e(2,jj)) * (sshn_e(3,jj) - hbdy_w(jj)) ) 415 END DO 416 ENDIF 417 ! 418 IF((nbondi == 1).OR.(nbondi == 2)) THEN 415 ua_e(ji,jj) = ubdy(ji,jj) * hur_e(ji,jj) 416 END DO 417 END DO 418 ! Characteristics method (only at boundary point): 419 ! istart = 2 420 ! iend = 2 421 ! DO ji = mi0(istart), mi1(iend) 422 ! DO jj=1,jpj 423 ! ua_e(ji,jj) = 0.5_wp * ( ubdy(ji,jj) * hur_e(ji,jj) + ua_e(ji+1,jj) & 424 ! & - sqrt(grav * hur_e(ji,jj)) * (sshn_e(ji+1,jj) - hbdy(ji,jj)) ) 425 ! END DO 426 ! END DO 427 ! 428 !--- East ---! 429 istart = jpiglo-nbghostcells 430 iend = jpiglo-1 431 DO ji = mi0(istart), mi1(iend) 419 432 DO jj=1,jpj 420 va_e(nlci-nbghostcells:nlci-1,jj) = vbdy_e(1:nbghostcells,jj) * hvr_e(nlci-nbghostcells:nlci-1,jj) 433 va_e(ji,jj) = vbdy(ji,jj) * hvr_e(ji,jj) 434 END DO 435 END DO 436 istart = jpiglo-nbghostcells-1 437 iend = jpiglo-2 438 DO ji = mi0(istart), mi1(iend) 439 DO jj=1,jpj 440 ua_e(ji,jj) = ubdy(ji,jj) * hur_e(ji,jj) 441 END DO 442 END DO 443 ! Characteristics method (only at boundary point): 444 ! istart = jpiglo-2 445 ! iend = jpiglo-2 446 ! DO ji = mi0(istart), mi1(iend) 447 ! DO jj=1,jpj 448 ! ua_e(ji,jj) = 0.5_wp * ( ubdy(ji,jj) * hur_e(ji,jj) + ua_e(ji-1,jj) & 449 ! & + sqrt(grav * hur_e(ji,jj)) * (sshn_e(ji,jj) - hbdy(ji+1,jj)) ) 450 ! END DO 451 ! END DO 452 ! 453 !--- South ---! 454 jstart = 2 455 jend = nbghostcells+1 456 DO jj = mj0(jstart), mj1(jend) 457 DO ji=1,jpi 458 ua_e(ji,jj) = ubdy(ji,jj) * hur_e(ji,jj) 421 459 ! Specified fluxes: 422 ua_e(nlci-nbghostcells-1:nlci-2,jj) = ubdy_e(1:nbghostcells,jj) * hur_e(nlci-nbghostcells-1:nlci-2,jj) 423 ! Characteristics method (only if ghostcells=1): 424 !alt ua_e(nlci-2,jj) = 0.5_wp * ( ubdy_e(jj) * hur_e(nlci-2,jj) + ua_e(nlci-3,jj) & 425 !alt & + sqrt(grav * hur_e(nlci-2,jj)) * (sshn_e(nlci-2,jj) - hbdy_e(jj)) ) 426 END DO 427 ENDIF 428 ! 429 IF((nbondj == -1).OR.(nbondj == 2)) THEN 460 va_e(ji,jj) = vbdy(ji,jj) * hvr_e(ji,jj) 461 END DO 462 END DO 463 ! Characteristics method (only at boundary point): 464 ! jstart = 2 465 ! jend = 2 466 ! DO jj = mj0(jstart), mj1(jend) 467 ! DO ji=1,jpi 468 ! va_e(ji,jj) = 0.5_wp * ( vbdy(ji,jj) * hvr_e(ji,jj) + va_e(ji,jj+1) & 469 ! & - sqrt(grav * hvr_e(ji,jj)) * (sshn_e(ji,jj+1) - hbdy(ji,jj)) ) 470 ! END DO 471 ! END DO 472 ! 473 !--- North ---! 474 jstart = jpjglo-nbghostcells 475 jend = jpjglo-1 476 DO jj = mj0(jstart), mj1(jend) 430 477 DO ji=1,jpi 431 ua_e(ji,2:nbghostcells+1) = ubdy_s(ji,1:nbghostcells) * hur_e(ji,2:nbghostcells+1) 432 ! Specified fluxes: 433 va_e(ji,2:nbghostcells+1) = vbdy_s(ji,1:nbghostcells) * hvr_e(ji,2:nbghostcells+1) 434 ! Characteristics method (only if ghostcells=1): 435 !alt va_e(ji,2) = 0.5_wp * ( vbdy_s(ji) * hvr_e(ji,2) + va_e(ji,3) & 436 !alt & - sqrt(grav * hvr_e(ji,2)) * (sshn_e(ji,3) - hbdy_s(ji)) ) 437 END DO 438 ENDIF 439 ! 440 IF((nbondj == 1).OR.(nbondj == 2)) THEN 478 ua_e(ji,jj) = ubdy(ji,jj) * hur_e(ji,jj) 479 END DO 480 END DO 481 jstart = jpjglo-nbghostcells-1 482 jend = jpjglo-2 483 DO jj = mj0(jstart), mj1(jend) 441 484 DO ji=1,jpi 442 ua_e(ji,nlcj-nbghostcells:nlcj-1) = ubdy_n(ji,1:nbghostcells) * hur_e(ji,nlcj-nbghostcells:nlcj-1) 443 ! Specified fluxes: 444 va_e(ji,nlcj-nbghostcells-1:nlcj-2) = vbdy_n(ji,1:nbghostcells) * hvr_e(ji,nlcj-nbghostcells-1:nlcj-2) 445 ! Characteristics method (only if ghostcells=1): 446 !alt va_e(ji,nlcj-2) = 0.5_wp * ( vbdy_n(ji) * hvr_e(ji,nlcj-2) + va_e(ji,nlcj-3) & 447 !alt & + sqrt(grav * hvr_e(ji,nlcj-2)) * (sshn_e(ji,nlcj-2) - hbdy_n(ji)) ) 448 END DO 449 ENDIF 485 va_e(ji,jj) = vbdy(ji,jj) * hvr_e(ji,jj) 486 END DO 487 END DO 488 ! Characteristics method (only at boundary point): 489 ! jstart = jpjglo-2 490 ! jend = jpjglo-2 491 ! DO jj = mj0(jstart), mj1(jend) 492 ! DO ji=1,jpi 493 ! va_e(ji,jj) = 0.5_wp * ( vbdy(ji,jj) * hvr_e(ji,jj) + va_e(ji,jj-1) & 494 ! & + sqrt(grav * hvr_e(ji,jj)) * (sshn_e(ji,jj) - hbdy(ji,jj+1)) ) 495 ! END DO 496 ! END DO 450 497 ! 451 498 END SUBROUTINE Agrif_dyn_ts … … 485 532 ELSE ! Linear interpolation 486 533 bdy_tinterp = 0 487 ubdy_w(:,:) = 0._wp ; vbdy_w(:,:) = 0._wp 488 ubdy_e(:,:) = 0._wp ; vbdy_e(:,:) = 0._wp 489 ubdy_n(:,:) = 0._wp ; vbdy_n(:,:) = 0._wp 490 ubdy_s(:,:) = 0._wp ; vbdy_s(:,:) = 0._wp 534 ubdy(:,:) = 0._wp ; vbdy(:,:) = 0._wp 491 535 CALL Agrif_Bc_variable( unb_id, procname=interpunb ) 492 536 CALL Agrif_Bc_variable( vnb_id, procname=interpvnb ) … … 503 547 INTEGER, INTENT(in) :: kt 504 548 ! 505 INTEGER :: ji, jj, indx, indy 549 INTEGER :: ji, jj 550 INTEGER :: istart, iend, jstart, jend 506 551 !!---------------------------------------------------------------------- 507 552 ! … … 516 561 ! 517 562 ! --- West --- ! 518 IF((nbondi == -1).OR.(nbondi == 2)) THEN 519 indx = 1+nbghostcells 563 istart = 2 564 iend = 1 + nbghostcells 565 DO ji = mi0(istart), mi1(iend) 520 566 DO jj = 1, jpj 521 DO ji = 2, indx 522 ssha(ji,jj) = hbdy_w(ji-1,jj) 523 ENDDO 567 ssha(ji,jj) = hbdy(ji,jj) 524 568 ENDDO 525 END IF569 ENDDO 526 570 ! 527 571 ! --- East --- ! 528 IF((nbondi == 1).OR.(nbondi == 2)) THEN 529 indx = nlci-nbghostcells 572 istart = jpiglo - nbghostcells 573 iend = jpiglo - 1 574 DO ji = mi0(istart), mi1(iend) 530 575 DO jj = 1, jpj 531 DO ji = indx, nlci-1 532 ssha(ji,jj) = hbdy_e(ji-indx+1,jj) 533 ENDDO 576 ssha(ji,jj) = hbdy(ji,jj) 534 577 ENDDO 535 END IF578 ENDDO 536 579 ! 537 580 ! --- South --- ! 538 IF((nbondj == -1).OR.(nbondj == 2)) THEN 539 indy = 1+nbghostcells 540 DO jj = 2, indy 541 DO ji = 1, jpi 542 ssha(ji,jj) = hbdy_s(ji,jj-1) 543 ENDDO 581 jstart = 2 582 jend = 1 + nbghostcells 583 DO jj = mj0(jstart), mj1(jend) 584 DO ji = 1, jpi 585 ssha(ji,jj) = hbdy(ji,jj) 544 586 ENDDO 545 END IF587 ENDDO 546 588 ! 547 589 ! --- North --- ! 548 IF((nbondj == 1).OR.(nbondj == 2)) THEN 549 indy = nlcj-nbghostcells 550 DO jj = indy, nlcj-1 551 DO ji = 1, jpi 552 ssha(ji,jj) = hbdy_n(ji,jj-indy+1) 553 ENDDO 590 jstart = jpjglo - nbghostcells 591 jend = jpjglo - 1 592 DO jj = mj0(jstart), mj1(jend) 593 DO ji = 1, jpi 594 ssha(ji,jj) = hbdy(ji,jj) 554 595 ENDDO 555 END IF596 ENDDO 556 597 ! 557 598 END SUBROUTINE Agrif_ssh … … 564 605 INTEGER, INTENT(in) :: jn 565 606 !! 566 INTEGER :: ji, jj , indx, indy567 !!----------------------------------------------------------------------568 !! clem ghost (starting at i,j=1 is important I think otherwise you introduce a grad(ssh)/=0 at point 2)607 INTEGER :: ji, jj 608 INTEGER :: istart, iend, jstart, jend 609 !!---------------------------------------------------------------------- 569 610 ! 570 611 IF( Agrif_Root() ) RETURN 571 612 ! 572 613 ! --- West --- ! 573 IF((nbondi == -1).OR.(nbondi == 2)) THEN 574 indx = 1+nbghostcells 614 istart = 2 615 iend = 1+nbghostcells 616 DO ji = mi0(istart), mi1(iend) 575 617 DO jj = 1, jpj 576 DO ji = 2, indx 577 ssha_e(ji,jj) = hbdy_w(ji-1,jj) 578 ENDDO 618 ssha_e(ji,jj) = hbdy(ji,jj) 579 619 ENDDO 580 END IF620 ENDDO 581 621 ! 582 622 ! --- East --- ! 583 IF((nbondi == 1).OR.(nbondi == 2)) THEN 584 indx = nlci-nbghostcells 623 istart = jpiglo - nbghostcells 624 iend = jpiglo - 1 625 DO ji = mi0(istart), mi1(iend) 585 626 DO jj = 1, jpj 586 DO ji = indx, nlci-1 587 ssha_e(ji,jj) = hbdy_e(ji-indx+1,jj) 588 ENDDO 627 ssha_e(ji,jj) = hbdy(ji,jj) 589 628 ENDDO 590 END IF629 ENDDO 591 630 ! 592 631 ! --- South --- ! 593 IF((nbondj == -1).OR.(nbondj == 2)) THEN 594 indy = 1+nbghostcells 595 DO jj = 2, indy 596 DO ji = 1, jpi 597 ssha_e(ji,jj) = hbdy_s(ji,jj-1) 598 ENDDO 632 jstart = 2 633 jend = 1+nbghostcells 634 DO jj = mj0(jstart), mj1(jend) 635 DO ji = 1, jpi 636 ssha_e(ji,jj) = hbdy(ji,jj) 599 637 ENDDO 600 END IF638 ENDDO 601 639 ! 602 640 ! --- North --- ! 603 IF((nbondj == 1).OR.(nbondj == 2)) THEN 604 indy = nlcj-nbghostcells 605 DO jj = indy, nlcj-1 606 DO ji = 1, jpi 607 ssha_e(ji,jj) = hbdy_n(ji,jj-indy+1) 608 ENDDO 641 jstart = jpjglo - nbghostcells 642 jend = jpjglo - 1 643 DO jj = mj0(jstart), mj1(jend) 644 DO ji = 1, jpi 645 ssha_e(ji,jj) = hbdy(ji,jj) 609 646 ENDDO 610 END IF647 ENDDO 611 648 ! 612 649 END SUBROUTINE Agrif_ssh_ts … … 848 885 END SUBROUTINE interptsn 849 886 850 SUBROUTINE interpsshn( ptab, i1, i2, j1, j2, before , nb, ndir)887 SUBROUTINE interpsshn( ptab, i1, i2, j1, j2, before ) 851 888 !!---------------------------------------------------------------------- 852 889 !! *** ROUTINE interpsshn *** … … 855 892 REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: ptab 856 893 LOGICAL , INTENT(in ) :: before 857 INTEGER , INTENT(in ) :: nb , ndir 858 ! 859 LOGICAL :: western_side, eastern_side,northern_side,southern_side 894 ! 860 895 !!---------------------------------------------------------------------- 861 896 ! … … 863 898 ptab(i1:i2,j1:j2) = sshn(i1:i2,j1:j2) 864 899 ELSE 865 western_side = (nb == 1).AND.(ndir == 1) 866 eastern_side = (nb == 1).AND.(ndir == 2) 867 southern_side = (nb == 2).AND.(ndir == 1) 868 northern_side = (nb == 2).AND.(ndir == 2) 869 !! clem ghost 870 IF(western_side) hbdy_w(1:nbghostcells,j1:j2) = ptab(i1:i2,j1:j2) * tmask(i1:i2,j1:j2,1) 871 IF(eastern_side) hbdy_e(1:nbghostcells,j1:j2) = ptab(i1:i2,j1:j2) * tmask(i1:i2,j1:j2,1) 872 IF(southern_side) hbdy_s(i1:i2,1:nbghostcells) = ptab(i1:i2,j1:j2) * tmask(i1:i2,j1:j2,1) 873 IF(northern_side) hbdy_n(i1:i2,1:nbghostcells) = ptab(i1:i2,j1:j2) * tmask(i1:i2,j1:j2,1) 900 hbdy(i1:i2,j1:j2) = ptab(i1:i2,j1:j2) * tmask(i1:i2,j1:j2,1) 874 901 ENDIF 875 902 ! … … 1045 1072 END SUBROUTINE interpvn 1046 1073 1047 SUBROUTINE interpunb( ptab, i1, i2, j1, j2, before , nb, ndir)1074 SUBROUTINE interpunb( ptab, i1, i2, j1, j2, before) 1048 1075 !!---------------------------------------------------------------------- 1049 1076 !! *** ROUTINE interpunb *** … … 1052 1079 REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: ptab 1053 1080 LOGICAL , INTENT(in ) :: before 1054 INTEGER , INTENT(in ) :: nb , ndir1055 1081 ! 1056 1082 INTEGER :: ji, jj 1057 1083 REAL(wp) :: zrhoy, zrhot, zt0, zt1, ztcoeff 1058 LOGICAL :: western_side, eastern_side,northern_side,southern_side1059 1084 !!---------------------------------------------------------------------- 1060 1085 ! … … 1062 1087 ptab(i1:i2,j1:j2) = e2u(i1:i2,j1:j2) * hu_n(i1:i2,j1:j2) * un_b(i1:i2,j1:j2) 1063 1088 ELSE 1064 western_side = (nb == 1).AND.(ndir == 1)1065 eastern_side = (nb == 1).AND.(ndir == 2)1066 southern_side = (nb == 2).AND.(ndir == 1)1067 northern_side = (nb == 2).AND.(ndir == 2)1068 1089 zrhoy = Agrif_Rhoy() 1069 1090 zrhot = Agrif_rhot() … … 1082 1103 ENDIF 1083 1104 ! 1084 IF(western_side) ubdy_w(1:nbghostcells,j1:j2) = ubdy_w(1:nbghostcells,j1:j2) + ztcoeff * ptab(i1:i2,j1:j2) 1085 IF(eastern_side) ubdy_e(1:nbghostcells,j1:j2) = ubdy_e(1:nbghostcells,j1:j2) + ztcoeff * ptab(i1:i2,j1:j2) 1086 IF(southern_side) ubdy_s(i1:i2,1:nbghostcells) = ubdy_s(i1:i2,1:nbghostcells) + ztcoeff * ptab(i1:i2,j1:j2) 1087 IF(northern_side) ubdy_n(i1:i2,1:nbghostcells) = ubdy_n(i1:i2,1:nbghostcells) + ztcoeff * ptab(i1:i2,j1:j2) 1105 ubdy(i1:i2,j1:j2) = ubdy(i1:i2,j1:j2) + ztcoeff * ptab(i1:i2,j1:j2) 1088 1106 ! 1089 1107 IF( bdy_tinterp == 0 .OR. bdy_tinterp == 2) THEN 1090 IF(western_side) ubdy_w(1:nbghostcells,j1:j2) = ubdy_w(1:nbghostcells,j1:j2) / (zrhoy*e2u(i1:i2,j1:j2)) * umask(i1:i2,j1:j2,1) 1091 IF(eastern_side) ubdy_e(1:nbghostcells,j1:j2) = ubdy_e(1:nbghostcells,j1:j2) / (zrhoy*e2u(i1:i2,j1:j2)) * umask(i1:i2,j1:j2,1) 1092 IF(southern_side) ubdy_s(i1:i2,1:nbghostcells) = ubdy_s(i1:i2,1:nbghostcells) / (zrhoy*e2u(i1:i2,j1:j2)) * umask(i1:i2,j1:j2,1) 1093 IF(northern_side) ubdy_n(i1:i2,1:nbghostcells) = ubdy_n(i1:i2,1:nbghostcells) / (zrhoy*e2u(i1:i2,j1:j2)) * umask(i1:i2,j1:j2,1) 1108 ubdy(i1:i2,j1:j2) = ubdy(i1:i2,j1:j2) / (zrhoy*e2u(i1:i2,j1:j2)) * umask(i1:i2,j1:j2,1) 1094 1109 ENDIF 1095 1110 ENDIF … … 1098 1113 1099 1114 1100 SUBROUTINE interpvnb( ptab, i1, i2, j1, j2, before , nb, ndir)1115 SUBROUTINE interpvnb( ptab, i1, i2, j1, j2, before ) 1101 1116 !!---------------------------------------------------------------------- 1102 1117 !! *** ROUTINE interpvnb *** … … 1105 1120 REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: ptab 1106 1121 LOGICAL , INTENT(in ) :: before 1107 INTEGER , INTENT(in ) :: nb , ndir1108 1122 ! 1109 1123 INTEGER :: ji,jj 1110 1124 REAL(wp) :: zrhox, zrhot, zt0, zt1, ztcoeff 1111 LOGICAL :: western_side, eastern_side,northern_side,southern_side1112 1125 !!---------------------------------------------------------------------- 1113 1126 ! … … 1115 1128 ptab(i1:i2,j1:j2) = e1v(i1:i2,j1:j2) * hv_n(i1:i2,j1:j2) * vn_b(i1:i2,j1:j2) 1116 1129 ELSE 1117 western_side = (nb == 1).AND.(ndir == 1)1118 eastern_side = (nb == 1).AND.(ndir == 2)1119 southern_side = (nb == 2).AND.(ndir == 1)1120 northern_side = (nb == 2).AND.(ndir == 2)1121 1130 zrhox = Agrif_Rhox() 1122 1131 zrhot = Agrif_rhot() … … 1133 1142 ztcoeff = 1 1134 1143 ENDIF 1135 !! clem ghost 1136 IF(western_side) vbdy_w(1:nbghostcells,j1:j2) = vbdy_w(1:nbghostcells,j1:j2) + ztcoeff * ptab(i1:i2,j1:j2) 1137 IF(eastern_side) vbdy_e(1:nbghostcells,j1:j2) = vbdy_e(1:nbghostcells,j1:j2) + ztcoeff * ptab(i1:i2,j1:j2) 1138 IF(southern_side) vbdy_s(i1:i2,1:nbghostcells) = vbdy_s(i1:i2,1:nbghostcells) + ztcoeff * ptab(i1:i2,j1:j2) 1139 IF(northern_side) vbdy_n(i1:i2,1:nbghostcells) = vbdy_n(i1:i2,1:nbghostcells) + ztcoeff * ptab(i1:i2,j1:j2) 1144 vbdy(i1:i2,j1:j2) = vbdy(i1:i2,j1:j2) + ztcoeff * ptab(i1:i2,j1:j2) 1140 1145 ! 1141 1146 IF( bdy_tinterp == 0 .OR. bdy_tinterp == 2) THEN 1142 IF(western_side) vbdy_w(1:nbghostcells,j1:j2) = vbdy_w(1:nbghostcells,j1:j2) / (zrhox*e1v(i1:i2,j1:j2)) * vmask(i1:i2,j1:j2,1) 1143 IF(eastern_side) vbdy_e(1:nbghostcells,j1:j2) = vbdy_e(1:nbghostcells,j1:j2) / (zrhox*e1v(i1:i2,j1:j2)) * vmask(i1:i2,j1:j2,1) 1144 IF(southern_side) vbdy_s(i1:i2,1:nbghostcells) = vbdy_s(i1:i2,1:nbghostcells) / (zrhox*e1v(i1:i2,j1:j2)) * vmask(i1:i2,j1:j2,1) 1145 IF(northern_side) vbdy_n(i1:i2,1:nbghostcells) = vbdy_n(i1:i2,1:nbghostcells) / (zrhox*e1v(i1:i2,j1:j2)) * vmask(i1:i2,j1:j2,1) 1146 ENDIF 1147 vbdy(i1:i2,j1:j2) = vbdy(i1:i2,j1:j2) / (zrhox*e1v(i1:i2,j1:j2)) * vmask(i1:i2,j1:j2,1) 1148 ENDIF 1147 1149 ENDIF 1148 1150 ! … … 1150 1152 1151 1153 1152 SUBROUTINE interpub2b( ptab, i1, i2, j1, j2, before , nb, ndir)1154 SUBROUTINE interpub2b( ptab, i1, i2, j1, j2, before ) 1153 1155 !!---------------------------------------------------------------------- 1154 1156 !! *** ROUTINE interpub2b *** … … 1157 1159 REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: ptab 1158 1160 LOGICAL , INTENT(in ) :: before 1159 INTEGER , INTENT(in ) :: nb , ndir1160 1161 ! 1161 1162 INTEGER :: ji,jj 1162 REAL(wp) :: zrhot, zt0, zt1,zat 1163 LOGICAL :: western_side, eastern_side,northern_side,southern_side 1163 REAL(wp) :: zrhot, zt0, zt1, zat 1164 1164 !!---------------------------------------------------------------------- 1165 1165 IF( before ) THEN … … 1170 1170 ENDIF 1171 1171 ELSE 1172 western_side = (nb == 1).AND.(ndir == 1)1173 eastern_side = (nb == 1).AND.(ndir == 2)1174 southern_side = (nb == 2).AND.(ndir == 1)1175 northern_side = (nb == 2).AND.(ndir == 2)1176 1172 zrhot = Agrif_rhot() 1177 1173 ! Time indexes bounds for integration … … 1181 1177 zat = zrhot * ( zt1**2._wp * (-2._wp*zt1 + 3._wp) & 1182 1178 & - zt0**2._wp * (-2._wp*zt0 + 3._wp) ) 1183 !! clem ghost 1184 IF(western_side ) ubdy_w(1:nbghostcells,j1:j2) = zat * ptab(i1:i2,j1:j2) 1185 IF(eastern_side ) ubdy_e(1:nbghostcells,j1:j2) = zat * ptab(i1:i2,j1:j2) 1186 IF(southern_side) ubdy_s(i1:i2,1:nbghostcells) = zat * ptab(i1:i2,j1:j2) 1187 IF(northern_side) ubdy_n(i1:i2,1:nbghostcells) = zat * ptab(i1:i2,j1:j2) 1179 ! 1180 ubdy(i1:i2,j1:j2) = zat * ptab(i1:i2,j1:j2) 1188 1181 ENDIF 1189 1182 ! … … 1191 1184 1192 1185 1193 SUBROUTINE interpvb2b( ptab, i1, i2, j1, j2, before , nb, ndir)1186 SUBROUTINE interpvb2b( ptab, i1, i2, j1, j2, before ) 1194 1187 !!---------------------------------------------------------------------- 1195 1188 !! *** ROUTINE interpvb2b *** … … 1198 1191 REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: ptab 1199 1192 LOGICAL , INTENT(in ) :: before 1200 INTEGER , INTENT(in ) :: nb , ndir1201 1193 ! 1202 1194 INTEGER :: ji,jj 1203 1195 REAL(wp) :: zrhot, zt0, zt1,zat 1204 LOGICAL :: western_side, eastern_side,northern_side,southern_side1205 1196 !!---------------------------------------------------------------------- 1206 1197 ! … … 1212 1203 ENDIF 1213 1204 ELSE 1214 western_side = (nb == 1).AND.(ndir == 1)1215 eastern_side = (nb == 1).AND.(ndir == 2)1216 southern_side = (nb == 2).AND.(ndir == 1)1217 northern_side = (nb == 2).AND.(ndir == 2)1218 1205 zrhot = Agrif_rhot() 1219 1206 ! Time indexes bounds for integration … … 1224 1211 & - zt0**2._wp * (-2._wp*zt0 + 3._wp) ) 1225 1212 ! 1226 IF(western_side ) vbdy_w(1:nbghostcells,j1:j2) = zat * ptab(i1:i2,j1:j2) 1227 IF(eastern_side ) vbdy_e(1:nbghostcells,j1:j2) = zat * ptab(i1:i2,j1:j2) 1228 IF(southern_side) vbdy_s(i1:i2,1:nbghostcells) = zat * ptab(i1:i2,j1:j2) 1229 IF(northern_side) vbdy_n(i1:i2,1:nbghostcells) = zat * ptab(i1:i2,j1:j2) 1213 vbdy(i1:i2,j1:j2) = zat * ptab(i1:i2,j1:j2) 1230 1214 ENDIF 1231 1215 ! -
NEMO/branches/2019/dev_r10973_AGRIF-01_jchanut_small_jpi_jpj/src/NST/agrif_oce_sponge.F90
r10425 r11205 93 93 !!---------------------------------------------------------------------- 94 94 INTEGER :: ji, jj, ind1, ind2 95 INTEGER :: ispongearea 96 REAL(wp) :: z1_ spongearea95 INTEGER :: ispongearea, jspongearea 96 REAL(wp) :: z1_ispongearea, z1_jspongearea 97 97 REAL(wp), DIMENSION(jpi,jpj) :: ztabramp 98 98 !!---------------------------------------------------------------------- … … 104 104 105 105 ispongearea = 1 + nn_sponge_len * Agrif_irhox() 106 z1_spongearea = 1._wp / REAL( ispongearea ) 106 z1_ispongearea = 1._wp / REAL( ispongearea ) 107 jspongearea = 1 + nn_sponge_len * Agrif_irhoy() 108 z1_jspongearea = 1._wp / REAL( jspongearea ) 107 109 108 110 ztabramp(:,:) = 0._wp 109 111 110 112 ! --- West --- ! 111 IF( (nbondi == -1) .OR. (nbondi == 2) ) THEN 112 ind1 = 1+nbghostcells 113 ind2 = 1+nbghostcells + ispongearea 113 ind1 = 1+nbghostcells 114 ind2 = 1+nbghostcells + ispongearea 115 DO ji = mi0(ind1), mi1(ind2) 116 DO jj = 1, jpj 117 ztabramp(ji,jj) = REAL( ind2 - mig(ji) ) * z1_ispongearea 118 END DO 119 END DO 120 121 ! --- East --- ! 122 ind1 = jpiglo - nbghostcells - ispongearea 123 ind2 = jpiglo - nbghostcells 124 DO ji = mi0(ind1), mi1(ind2) 114 125 DO jj = 1, jpj 115 DO ji = ind1, ind2 116 ztabramp(ji,jj) = REAL( ind2 - ji ) * z1_spongearea * umask(ind1,jj,1) 117 END DO 126 ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL( mig(ji) - ind1 ) * z1_ispongearea) 118 127 ENDDO 119 ENDIF 120 121 ! --- East --- ! 122 IF( (nbondi == 1) .OR. (nbondi == 2) ) THEN 123 ind1 = nlci - nbghostcells - ispongearea 124 ind2 = nlci - nbghostcells 125 DO jj = 1, jpj 126 DO ji = ind1, ind2 127 ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL( ji - ind1 ) * z1_spongearea * umask(ind2-1,jj,1) ) 128 ENDDO 129 ENDDO 130 ENDIF 128 END DO 131 129 132 130 ! --- South --- ! 133 IF( (nbondj == -1) .OR. (nbondj == 2) ) THEN 134 ind1 = 1+nbghostcells 135 ind2 = 1+nbghostcells + ispongearea 136 DO jj = ind1, ind2 137 DO ji = 1, jpi 138 ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL( ind2 - jj ) * z1_spongearea * vmask(ji,ind1,1) ) 139 END DO 140 ENDDO 141 ENDIF 131 ind1 = 1+nbghostcells 132 ind2 = 1+nbghostcells + jspongearea 133 DO jj = mj0(ind1), mj1(ind2) 134 DO ji = 1, jpi 135 ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL( ind2 - mjg(jj) ) * z1_jspongearea) 136 END DO 137 END DO 142 138 143 139 ! --- North --- ! 144 IF( (nbondj == 1) .OR. (nbondj == 2) ) THEN 145 ind1 = nlcj - nbghostcells - ispongearea 146 ind2 = nlcj - nbghostcells 147 DO jj = ind1, ind2 148 DO ji = 1, jpi 149 ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL( jj - ind1 ) * z1_spongearea * vmask(ji,ind2-1,1) ) 150 END DO 151 ENDDO 152 ENDIF 140 ind1 = jpjglo - nbghostcells - jspongearea 141 ind2 = jpjglo - nbghostcells 142 DO jj = mj0(ind1), mj1(ind2) 143 DO ji = 1, jpi 144 ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL( mjg(jj) - ind1 ) * z1_jspongearea) 145 END DO 146 END DO 153 147 154 148 ENDIF -
NEMO/branches/2019/dev_r10973_AGRIF-01_jchanut_small_jpi_jpj/src/NST/agrif_user.F90
r10425 r11205 190 190 Agrif_UseSpecialValue = .TRUE. 191 191 CALL Agrif_Bc_variable(sshn_id,calledweight=1., procname=interpsshn ) 192 hbdy _w(:,:) = 0.e0 ; hbdy_e(:,:) = 0.e0 ; hbdy_n(:,:) = 0.e0 ; hbdy_s(:,:) = 0.e0192 hbdy(:,:) = 0.e0 193 193 ssha(:,:) = 0.e0 194 194 … … 199 199 CALL Agrif_Bc_variable(ub2b_interp_id,calledweight=1.,procname=interpub2b) 200 200 CALL Agrif_Bc_variable(vb2b_interp_id,calledweight=1.,procname=interpvb2b) 201 ubdy_w(:,:) = 0.e0 ; vbdy_w(:,:) = 0.e0 202 ubdy_e(:,:) = 0.e0 ; vbdy_e(:,:) = 0.e0 203 ubdy_n(:,:) = 0.e0 ; vbdy_n(:,:) = 0.e0 204 ubdy_s(:,:) = 0.e0 ; vbdy_s(:,:) = 0.e0 201 ubdy(:,:) = 0.e0 ; vbdy(:,:) = 0.e0 205 202 ENDIF 206 203 … … 736 733 ! 737 734 ! Check sponge length: 738 iminspon = MIN(FLOOR(REAL(jpiglo-4)/REAL(2*Agrif_irhox())), FLOOR(REAL(jpjglo-4)/REAL(2*Agrif_irhox())) )739 IF (lk_mpp) iminspon = MIN(iminspon,FLOOR(REAL(jpi-2)/REAL(Agrif_irhox())), FLOOR(REAL(jpj-2)/REAL(Agrif_irhox())) )740 IF (nn_sponge_len > iminspon) CALL ctl_stop('agrif sponge length is too large')735 ! iminspon = MIN(FLOOR(REAL(jpiglo-4)/REAL(2*Agrif_irhox())), FLOOR(REAL(jpjglo-4)/REAL(2*Agrif_irhox())) ) 736 ! IF (lk_mpp) iminspon = MIN(iminspon,FLOOR(REAL(jpi-2)/REAL(Agrif_irhox())), FLOOR(REAL(jpj-2)/REAL(Agrif_irhox())) ) 737 ! IF (nn_sponge_len > iminspon) CALL ctl_stop('agrif sponge length is too large') 741 738 ! 742 739 IF( agrif_oce_alloc() > 0 ) CALL ctl_warn('agrif agrif_oce_alloc: allocation of arrays failed') -
NEMO/branches/2019/dev_r10973_AGRIF-01_jchanut_small_jpi_jpj/src/OCE/DYN/dynspg_ts.F90
r10742 r11205 797 797 IF( ln_bdy .AND. ln_vol ) CALL bdy_vol2d( kt, jn, ua_e, va_e, zhup2_e, zhvp2_e ) 798 798 ! 799 #if defined key_agrif 800 ! Set fluxes during predictor step to ensure volume conservation 801 IF( .NOT.Agrif_Root() .AND. ln_bt_fw ) CALL agrif_dyn_ts( jn ) 802 #endif 799 803 zwx(:,:) = e2u(:,:) * ua_e(:,:) * zhup2_e(:,:) ! fluxes at jn+0.5 800 804 zwy(:,:) = e1v(:,:) * va_e(:,:) * zhvp2_e(:,:) 801 805 ! 802 #if defined key_agrif803 ! Set fluxes during predictor step to ensure volume conservation804 IF( .NOT.Agrif_Root() .AND. ln_bt_fw ) THEN805 IF((nbondi == -1).OR.(nbondi == 2)) THEN806 DO jj = 1, jpj807 zwx(2:nbghostcells+1,jj) = ubdy_w(1:nbghostcells,jj) * e2u(2:nbghostcells+1,jj)808 zwy(2:nbghostcells+1,jj) = vbdy_w(1:nbghostcells,jj) * e1v(2:nbghostcells+1,jj)809 END DO810 ENDIF811 IF((nbondi == 1).OR.(nbondi == 2)) THEN812 DO jj=1,jpj813 zwx(nlci-nbghostcells-1:nlci-2,jj) = ubdy_e(1:nbghostcells,jj) * e2u(nlci-nbghostcells-1:nlci-2,jj)814 zwy(nlci-nbghostcells :nlci-1,jj) = vbdy_e(1:nbghostcells,jj) * e1v(nlci-nbghostcells :nlci-1,jj)815 END DO816 ENDIF817 IF((nbondj == -1).OR.(nbondj == 2)) THEN818 DO ji=1,jpi819 zwy(ji,2:nbghostcells+1) = vbdy_s(ji,1:nbghostcells) * e1v(ji,2:nbghostcells+1)820 zwx(ji,2:nbghostcells+1) = ubdy_s(ji,1:nbghostcells) * e2u(ji,2:nbghostcells+1)821 END DO822 ENDIF823 IF((nbondj == 1).OR.(nbondj == 2)) THEN824 DO ji=1,jpi825 zwy(ji,nlcj-nbghostcells-1:nlcj-2) = vbdy_n(ji,1:nbghostcells) * e1v(ji,nlcj-nbghostcells-1:nlcj-2)826 zwx(ji,nlcj-nbghostcells :nlcj-1) = ubdy_n(ji,1:nbghostcells) * e2u(ji,nlcj-nbghostcells :nlcj-1)827 END DO828 ENDIF829 ENDIF830 #endif831 806 IF( ln_wd_il ) CALL wad_lmt_bt(zwx, zwy, sshn_e, zssh_frc, rdtbt) 832 807
Note: See TracChangeset
for help on using the changeset viewer.