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 26 for trunk/NEMO/OPA_SRC/OBC – NEMO

Changeset 26 for trunk/NEMO/OPA_SRC/OBC


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

CT : BUGFIX010 : Running problem is solved

File:
1 edited

Legend:

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

    r3 r26  
    1919   USE daymod          ! calendar 
    2020   USE in_out_manager  ! I/O logical units 
     21   USE lib_mpp         ! distribued memory computing 
    2122 
    2223 
     
    4243   SUBROUTINE obc_dta_uvt ( kt ) 
    4344      !!--------------------------------------------------------------------------- 
    44       !!                         SUBROUTINE obc_dta_uvt  
    45       !!                        ************************ 
    46       !! ** Purpose : 
    47       !!      Find the climatological  boundary arrays for the specified date,  
     45      !!                      ***  SUBROUTINE obc_dta_uvt  *** 
     46      !!                     
     47      !! ** Purpose :   Find the climatological  boundary arrays for the specified date,  
    4848      !!      Originally this routine interpolated between monthly fields 
    4949      !!      of a climatology. 
     
    5151      !!      and do not need to interpolate. 
    5252      !! 
    53       !! ** Method : 
    54       !!      Determine the current month from kdat, and interpolate for the 
    55       !!      current day. 
     53      !! ** Method  :   Determine the current month from kdat, and interpolate for 
     54      !!      the current day. 
    5655      !! 
    5756      !! History : 
     
    149148                        sedta(ij,jk,1) = sn(ji,jj,jk)*tmask(ji,jj,jk) 
    150149                        tedta(ij,jk,1) = tn(ji,jj,jk)*tmask(ji,jj,jk) 
    151                         uedta(ij,jk,1) = 0.1*umask(ji,jj,jk) 
     150                        uedta(ij,jk,1) = un(ji,jj,jk)*umask(ji,jj,jk) 
    152151                     END DO 
    153152                  END DO 
     
    240239            IF( nobc_dta == 0 )   THEN                ! initial state used 
    241240               !                                      ! ================== ! 
    242             DO ji = fs_niw0, fs_niw1 ! Vector opt. 
     241               DO ji = fs_niw0, fs_niw1 ! Vector opt. 
     242                  DO jk = 1, jpkm1 
     243                     DO jj = 1, jpj 
     244                        ij = jj -1 + njmpp 
     245                        swdta(ij,jk,1) = sn(ji,jj,jk)*tmask(ji,jj,jk) 
     246                        twdta(ij,jk,1) = tn(ji,jj,jk)*tmask(ji,jj,jk) 
     247                        uwdta(ij,jk,1) = un(ji,jj,jk)*umask(ji,jj,jk) 
     248                     END DO 
     249                  END DO 
     250               END DO 
     251    
    243252               DO jk = 1, jpkm1 
    244253                  DO jj = 1, jpj 
    245254                     ij = jj -1 + njmpp 
    246                      swdta(ij,jk,1) = sn(ji,jj,jk)*tmask(ji,jj,jk) 
    247                      twdta(ij,jk,1) = tn(ji,jj,jk)*tmask(ji,jj,jk) 
    248                      uwdta(ij,jk,1) = 0.1*umask(ji,jj,jk) 
     255                     sfow(jj,jk) = swdta(ij,jk,1)*twmsk(jj,jk) 
     256                     tfow(jj,jk) = twdta(ij,jk,1)*twmsk(jj,jk) 
     257                     ufow(jj,jk) = uwdta(ij,jk,1)*uwmsk(jj,jk) 
    249258                  END DO 
    250259               END DO 
    251             END DO 
    252  
    253             DO jk = 1, jpkm1 
    254                DO jj = 1, jpj 
    255                   ij = jj -1 + njmpp 
    256                   sfow(jj,jk) = swdta(ij,jk,1)*twmsk(jj,jk) 
    257                   tfow(jj,jk) = twdta(ij,jk,1)*twmsk(jj,jk) 
    258                   ufow(jj,jk) = uwdta(ij,jk,1)*uwmsk(jj,jk) 
    259                END DO 
    260             END DO 
    261260               !                                      ! =================== ! 
    262261            ELSE                                      ! read in obceast.dta 
    263262               !                                      ! =================== ! 
    264             OPEN(UNIT   = inum,       & 
     263               OPEN(UNIT   = inum,      & 
    265264                 IOSTAT = ios,          & 
    266265                 FILE   ='obcwest.dta', & 
     
    268267                 ACCESS ='DIRECT',      & 
    269268                 RECL   = 4096 ) 
    270             IF( ios > 0 ) THEN 
    271                IF(lwp) WRITE(numout,*) 'obc_dta_uvt: Pbm to OPEN the obcwest.dta file ' 
    272                IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' 
    273                nstop = nstop + 1 
    274             END IF 
    275             READ(inum,REC=1) clversion, clcom,irecl 
    276             CLOSE(inum) 
    277             IF(lwp) WRITE(numout,*)'       ' 
    278             IF(lwp) WRITE(numout,*)'         opening obcwest.dta  with irecl=',irecl 
    279             OPEN(UNIT   = inum,       & 
     269               IF( ios > 0 ) THEN 
     270                  IF(lwp) WRITE(numout,*) 'obc_dta_uvt: Pbm to OPEN the obcwest.dta file ' 
     271                  IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' 
     272                  nstop = nstop + 1 
     273               END IF 
     274               READ(inum,REC=1) clversion, clcom,irecl 
     275               CLOSE(inum) 
     276               IF(lwp) WRITE(numout,*)'       ' 
     277               IF(lwp) WRITE(numout,*)'         opening obcwest.dta  with irecl=',irecl 
     278               OPEN(UNIT   = inum,      & 
    280279                 IOSTAT = ios,          & 
    281280                 FILE   ='obcwest.dta', & 
     
    283282                 ACCESS ='DIRECT',      & 
    284283                 RECL   = irecl ) 
    285             IF( ios > 0 ) THEN 
    286                IF(lwp) WRITE(numout,*) 'obc_dta_uvt: Pbm to OPEN the obcwest.dta file ' 
    287                IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' 
    288                nstop = nstop + 1 
    289             END IF 
    290  
    291             ! ... Read datafile and set temperature, salinity and normal velocity 
    292             ! ... initialise the swdta, twdta arrays 
    293             ! ... index 1 refer to before, 2 to after 
    294             DO jk = 1, jpkm1 
    295                irec = 2 +  (jk -1)* jpf 
    296                READ(inum,REC=irec  )((swdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 
    297                READ(inum,REC=irec+1)((twdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 
    298                READ(inum,REC=irec+2)((uwdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 
    299                DO jj = 1, jpj 
    300                   ij = jj -1 + njmpp 
    301                   sfow(jj,jk) = swdta(ij,jk,1)*twmsk(jj,jk) 
    302                   tfow(jj,jk) = twdta(ij,jk,1)*twmsk(jj,jk) 
    303                END DO 
    304             END DO 
    305             CLOSE(inum) 
     284               IF( ios > 0 ) THEN 
     285                  IF(lwp) WRITE(numout,*) 'obc_dta_uvt: Pbm to OPEN the obcwest.dta file ' 
     286                  IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' 
     287                  nstop = nstop + 1 
     288               END IF 
     289 
     290               ! ... Read datafile and set temperature, salinity and normal velocity 
     291               ! ... initialise the swdta, twdta arrays 
     292               ! ... index 1 refer to before, 2 to after 
     293               DO jk = 1, jpkm1 
     294                  irec = 2 +  (jk -1)* jpf 
     295                  READ(inum,REC=irec  )((swdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 
     296                  READ(inum,REC=irec+1)((twdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 
     297                  READ(inum,REC=irec+2)((uwdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 
     298                  DO jj = 1, jpj 
     299                     ij = jj -1 + njmpp 
     300                     sfow(jj,jk) = swdta(ij,jk,1)*twmsk(jj,jk) 
     301                     tfow(jj,jk) = twdta(ij,jk,1)*twmsk(jj,jk) 
     302                  END DO 
     303               END DO 
     304               CLOSE(inum) 
    306305 
    307306#if ! defined key_dynspg_fsc 
    308             ! ... Rigid lid case: make sure uwdta is baroclinic velocity 
    309             ! ... In rigid lid case uwdta needs to be the baroclinic component. 
    310  
    311             CALL obc_cli( uwdta, ucliw, fs_niw0, fs_niw1, 0, jpj, njmpp )   
     307               ! ... Rigid lid case: make sure uwdta is baroclinic velocity 
     308               ! ... In rigid lid case uwdta needs to be the baroclinic component. 
     309 
     310               CALL obc_cli( uwdta, ucliw, fs_niw0, fs_niw1, 0, jpj, njmpp )   
    312311 
    313312# endif 
    314             ! ... Set normal velocity (on niw0, niw1 <=> jpiwob) 
    315             DO jk = 1, jpkm1 
    316                 DO jj = 1, jpj 
    317                    ij = jj -1 + njmpp 
    318                    ufow(jj,jk) = uwdta(ij,jk,1)*uwmsk(jj,jk) 
    319                 END DO 
    320             END DO 
     313               ! ... Set normal velocity (on niw0, niw1 <=> jpiwob) 
     314               DO jk = 1, jpkm1 
     315                   DO jj = 1, jpj 
     316                      ij = jj -1 + njmpp 
     317                      ufow(jj,jk) = uwdta(ij,jk,1)*uwmsk(jj,jk) 
     318                   END DO 
     319               END DO 
    321320            ENDIF 
    322321         ENDIF 
     
    332331            vndta(:,:,1) = 0.e0 
    333332 
    334             OPEN(UNIT   = inum,        & 
     333            OPEN(UNIT   = inum,          & 
    335334                 IOSTAT = ios,           & 
    336335                 FILE   ='obcnorth.dta', & 
     
    528527      !! * Arguments 
    529528      INTEGER,INTENT(in) :: kt  
    530       WRITE(*,*) kt 
     529      WRITE(*,*) 'obc_dta_psi: You should not have seen this print! error?', kt 
    531530   END SUBROUTINE obc_dta_psi 
    532531#else 
     
    567566      !! * Local declarations 
    568567      INTEGER ::   ji, jj, jnic, jip         ! dummy loop indices 
     568      INTEGER ::   inum = 11                 ! temporary logical unit 
    569569      INTEGER ::   ip, ii, ij, iii, ijj 
    570570      INTEGER ::   kbsfstart 
     
    622622         END DO 
    623623      END IF 
    624 # if defined key_mpp 
    625       CALL mpprisl( gcbic, 3 ) 
    626 # endif 
     624 
     625      IF( lk_mpp )   CALL mpp_isl( gcbic, 3 ) 
    627626 
    628627      ! 3. Update the climatological barotropic function at the boundary  
     
    711710   SUBROUTINE obc_dta_uvt( kt )             ! Empty routine 
    712711      INTEGER, INTENT (in) :: kt 
    713       WRITE(*,*) kt 
     712      WRITE(*,*) 'obc_dta_uvt: You should not have seen this print! error?', kt 
    714713   END SUBROUTINE obc_dta_uvt 
    715714 
    716715   SUBROUTINE obc_dta_psi( kt )             ! Empty routine 
    717716      INTEGER, INTENT (in) :: kt 
    718       WRITE(*,*) kt 
     717      WRITE(*,*) 'obc_dta_psi: You should not have seen this print! error?', kt 
    719718   END SUBROUTINE obc_dta_psi 
    720719 
Note: See TracChangeset for help on using the changeset viewer.