New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 32 for trunk/NEMO/OPA_SRC/OBC/obcspg.F90 – NEMO

Ignore:
Timestamp:
2004-02-17T10:20:15+01:00 (21 years ago)
Author:
opalod
Message:

CT : UPDATE001 : First major NEMO update

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/OBC/obcspg.F90

    r3 r32  
    55   !!                      open boundary 
    66   !!====================================================================== 
    7 #if defined key_obc && defined key_dynspg_rl 
     7#if   defined key_obc   &&  defined key_dynspg_rl 
    88   !!---------------------------------------------------------------------- 
    99   !!   'key_obc'    and                            Open Boundary Condition 
     
    8686      !!---------------------------------------------------------------------- 
    8787 
    88       ! 0. Local constant initialization 
    89       ! -------------------------------- 
    90  
    91       IF( kt == nit000 .OR. ln_rstart ) THEN 
     88      IF( kt == nit000 .OR. ln_rstart ) THEN      ! Initialization 
    9289         ! ... Boundary restoring coefficient 
    9390         rtaue = 2. * rdt / rdpeob 
     
    10097         rtaunin = 2. * rdt / rdpnin 
    10198         rtausin = 2. * rdt / rdpsin  
    102       END IF 
    103  
    104       ! ... right hand side of the barotropic elliptic equation 
     99      ENDIF 
     100 
     101      ! right hand side of the barotropic elliptic equation 
     102      ! --------------------------------------------------- 
     103 
     104      ! Isolated coastline contribution to the RHS of the barotropic Eq. 
    105105      gcbob(:,:) = 0.e0 
    106  
    107       ! 1. Isolated coastline contribution to the RHS of the barotropic Eq. 
    108       ! ------------------------------------------------------------------- 
    109106      DO jnic = 1, nbobc-1 
    110          DO jj = 1, jpj 
    111             DO ji = 1, jpi 
    112                gcbob(ji,jj) = gcbob(ji,jj) + gcfobc(ji,jj,jnic) * gcbic(jnic) 
    113             END DO 
    114          END DO 
     107         gcbob(:,:) = gcbob(:,:) + gcfobc(:,:,jnic) * gcbic(jnic) 
    115108      END DO 
    116109 
    117       ! 2. East open boundary 
    118       ! --------------------- 
    119  
    120       IF( lpeastobc ) THEN 
    121          CALL obc_spg_east( kt ) 
    122       END IF 
    123  
    124       ! 3. West open boundary 
    125       ! --------------------- 
    126  
    127       IF( lpwestobc ) THEN 
    128          CALL obc_spg_west( kt ) 
    129       END IF 
    130          
    131       ! 4. North open boundary 
    132       ! ---------------------- 
    133  
    134       IF( lpnorthobc ) THEN 
    135          CALL obc_spg_north( kt ) 
    136       END IF 
    137   
    138       ! 5. South open boundary 
    139       ! ---------------------- 
    140  
    141       IF( lpsouthobc ) THEN 
    142          CALL obc_spg_south( kt ) 
    143       END IF 
    144   
    145 # if defined key_mpp 
    146       CALL mpp_lnk_2d( gcbob, 'G', 1. ) 
    147 # endif 
     110      IF( lpeastobc  )   CALL obc_spg_east ( kt )    ! East open boundary 
     111 
     112      IF( lpwestobc  )   CALL obc_spg_west ( kt )    ! West open boundary 
     113 
     114      IF( lpnorthobc )   CALL obc_spg_north( kt )    ! North open boundary 
     115 
     116      IF( lpsouthobc )   CALL obc_spg_south( kt )    ! South open boundary 
     117 
     118      IF( lk_mpp )   CALL lbc_lnk( gcbob, 'G', 1. ) 
    148119  
    149120   END SUBROUTINE obc_spg 
    150121 
     122 
    151123   SUBROUTINE obc_spg_east ( kt ) 
    152124      !!------------------------------------------------------------------------------ 
    153       !!                     SUBROUTINE obc_spg_east 
    154       !!                    ************************* 
    155       !! ** Purpose : 
    156       !!      Apply the radiation algorithm on east OBC stream function. 
    157       !!      If the logical lfbceast is .TRUE., there is no radiation but only fixed OBC 
     125      !!                ***  SUBROUTINE obc_spg_east  *** 
     126      !!                  
     127      !! ** Purpose :   Apply the radiation algorithm on east OBC stream function. 
     128      !!      If lfbceast=T , there is no radiation but only fixed OBC 
    158129      !! 
    159130      !!  History : 
     
    169140      !! * Local declarations 
    170141      INTEGER ::   ij 
    171  
    172142      REAL(wp) ::   z2dtr, ztau, zin 
    173143      REAL(wp) ::   z05cx, zdt, z4nor2, z2dx, z2dy 
    174  
    175       !!------------------------------------------------------------------------------ 
    176       !!  OPA 8.5, LODYC-IPSL (2002) 
    177144      !!------------------------------------------------------------------------------ 
    178145 
     
    229196                  IF(lwp) WRITE(numout,*)' PB dans obc_spg_east au pt ',jj,' : z4nor=0' 
    230197                  z4nor2 = 0.001 
    231                END IF 
     198               ENDIF 
    232199               z05cx = zdt * z2dx / z4nor2 * bmask(ji,jj) 
    233200               z05cx = z05cx / e1v(ji+1,jj) 
     
    249216         END DO 
    250217 
    251       END IF 
    252 # if defined key_mpp 
    253       CALL mppobc(bsfeob,jpjed,jpjef,jpieob-1,1,2,jpj) 
    254 # endif 
     218      ENDIF 
     219      IF( lk_mpp )   CALL mppobc(bsfeob,jpjed,jpjef,jpieob-1,1,2,jpj) 
     220 
    255221 
    256222      ! 3. right hand side of the barotropic elliptic equation 
     
    258224  
    259225      IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN 
    260          z2dtr=1./rdt 
     226         z2dtr = 1.0 / rdt 
    261227      ELSE 
    262          z2dtr=1./2./rdt 
    263       END IF 
     228         z2dtr = 0.5 / rdt 
     229      ENDIF 
    264230      DO ji = fs_nie0-1, fs_nie1-1 ! Vector opt. 
    265231         DO jj = nje0m1, nje1  
     
    351317                  IF(lwp) WRITE(numout,*)' PB dans obc_spg_west au pt ',jj,' : z4nor =0' 
    352318                  z4nor2=0.0001 
    353                END IF 
     319               ENDIF 
    354320               z05cx = zdt * z2dx / z4nor2 * bmask(ji,jj) 
    355321               z05cx = z05cx / e1v(ji,jj) 
     
    368334         END DO 
    369335 
    370       END IF 
    371 # if defined key_mpp  
    372       CALL mppobc(bsfwob,jpjwd,jpjwf,jpiwob+1,1,2,jpj)  
    373 # endif  
     336      ENDIF 
     337      IF( lk_mpp )   CALL mppobc(bsfwob,jpjwd,jpjwf,jpiwob+1,1,2,jpj)  
     338 
    374339 
    375340      ! 3. right hand side of the barotropic elliptic equation 
     
    377342 
    378343      IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN 
    379          z2dtr=1./rdt 
     344         z2dtr = 1.0 / rdt 
    380345      ELSE 
    381          z2dtr=1./2./rdt 
    382       END IF 
     346         z2dtr = 0.5 / rdt 
     347      ENDIF 
    383348      DO ji = fs_niw0+1, fs_niw1+1 ! Vector opt. 
    384349         DO jj = njw0m1, njw1 
     
    392357   SUBROUTINE obc_spg_north ( kt ) 
    393358      !!------------------------------------------------------------------------------ 
    394       !!                     SUBROUTINE obc_spg_north 
    395       !!                    ************************* 
    396       !! ** Purpose : 
    397       !!      Apply the radiation algorithm on north OBC stream function. 
    398       !!      If the logical lfbcnorth is .TRUE., there is no radiation but only fixed OBC 
     359      !!                 ***  SUBROUTINE obc_spg_north  *** 
     360      !!  
     361      !! ** Purpose :   Apply the radiation algorithm on north OBC stream function. 
     362      !!      If lfbcnorth=T, there is no radiation but only fixed OBC 
    399363      !! 
    400364      !!  History : 
     
    410374      !! * Local declarations 
    411375      INTEGER ::   ii 
    412  
    413376      REAL(wp) ::   z2dtr, ztau, zin 
    414377      REAL(wp) ::   z05cx, zdt, z4nor2, z2dx, z2dy 
    415  
    416       !!------------------------------------------------------------------------------ 
    417       !!  OPA 8.5, LODYC-IPSL (2002) 
    418378      !!------------------------------------------------------------------------------ 
    419379 
     
    475435               IF( z4nor2 == 0 ) THEN 
    476436                  IF(lwp) WRITE(numout,*)' PB dans obc_spg_north au pt',ji,' : z4nor =0' 
    477                END IF 
     437               ENDIF 
    478438               z05cx = zdt * z2dx / z4nor2 * bmask(ji,jj) 
    479439               z05cx = z05cx / e2u(ji,jj+1) 
     
    492452         END DO 
    493453 
    494       END IF 
    495 # if defined key_mpp 
    496       call mppobc(bsfnob,jpind,jpinf,jpjnob-1,1,1,jpi) 
    497 # endif 
     454      ENDIF 
     455      IF( lk_mpp )   CALL mppobc(bsfnob,jpind,jpinf,jpjnob-1,1,1,jpi) 
     456 
    498457 
    499458      ! 3. right hand side of the barotropic elliptic equation 
     
    501460 
    502461      IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN 
    503          z2dtr=1./rdt 
     462         z2dtr = 1.0 / rdt 
    504463      ELSE 
    505          z2dtr=1./2./rdt 
    506       END IF 
     464         z2dtr = 0.5 / rdt 
     465      ENDIF 
    507466      DO jj = fs_njn0-1, fs_njn1-1 ! Vector opt. 
    508467         DO ji = nin0m1, nin1 
     
    514473   END SUBROUTINE obc_spg_north 
    515474 
     475 
    516476   SUBROUTINE obc_spg_south ( kt ) 
    517477      !!------------------------------------------------------------------------------ 
    518       !!                     SUBROUTINE obc_spg_south 
    519       !!                    ************************* 
    520       !! ** Purpose : 
    521       !!      Apply the radiation algorithm on south OBC stream function. 
    522       !!      If the logical lfbcsouth is .TRUE., there is no radiation but only fixed OBC 
     478      !!                  ***  SUBROUTINE obc_spg_south  *** 
     479      !!                 
     480      !! ** Purpose :   Apply the radiation algorithm on south OBC stream function. 
     481      !!      If lfbcsouth=T, there is no radiation but only fixed OBC 
    523482      !! 
    524483      !!  History : 
     
    596555               IF( z4nor2 == 0 ) THEN 
    597556                  IF(lwp) WRITE(numout,*)' PB dans obc_spg_south au pt ',ji,' : z4nor =0' 
    598                END IF 
     557               ENDIF 
    599558               z05cx = zdt * z2dx / z4nor2 * bmask(ji,jj) 
    600559               z05cx = z05cx / e2u(ji,jj) 
     
    613572         END DO 
    614573 
    615       END IF 
    616 # if defined key_mpp 
    617       CALL mppobc(bsfsob,jpisd,jpisf,jpjsob+1,1,1,jpi) 
    618 # endif 
     574      ENDIF 
     575      IF( lk_mpp )   CALL mppobc(bsfsob,jpisd,jpisf,jpjsob+1,1,1,jpi) 
     576 
    619577  
    620578      ! 3. right hand side of the barotropic elliptic equation 
    621579      ! ------------------------------------------------------- 
    622580 
    623       IF( ( neuler == 0 ) .and. ( kt == nit000 ) ) THEN 
    624          z2dtr=1./rdt 
     581      IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN 
     582         z2dtr = 1.0 / rdt 
    625583      ELSE 
    626          z2dtr=1./2./rdt 
    627       END IF 
     584         z2dtr = 0.5 / rdt 
     585      ENDIF 
    628586      DO jj = fs_njs0+1, fs_njs1+1 ! Vector opt. 
    629587         DO ji = nis0m1, nis1  
     
    642600   SUBROUTINE obc_spg( kt )        ! Empty routine 
    643601      INTEGER, INTENT( in ) :: kt 
    644       WRITE(*,*) kt 
     602      WRITE(*,*) 'obc_spg: You should not have seen this print! error?', kt 
    645603   END SUBROUTINE obc_spg 
    646604#endif 
Note: See TracChangeset for help on using the changeset viewer.