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 5639 – NEMO

Changeset 5639


Ignore:
Timestamp:
2015-07-28T14:20:28+02:00 (9 years ago)
Author:
jamesharle
Message:

bug fix to fldread to allow code to run with bt vels properly

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/SBC/fldread.F90

    r5625 r5639  
    759759 
    760760      SELECT CASE( ipk ) 
    761       CASE(1)        ;   CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1    ), nrec ) 
    762       CASE DEFAULT   ;    
    763    
     761      CASE(1)        ;    
     762      CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1    ), nrec ) 
     763         IF (ipj==1) THEN ! we assume that this is an un-structured open boundary file 
     764            DO ib = 1, ipi 
     765              DO ik = 1, ipk 
     766                dta(ib,1,ik) =  dta_read(map%ptr(ib),1,ik) 
     767              END DO 
     768            END DO 
     769         ELSE ! we assume that this is a structured open boundary file 
     770            DO ib = 1, ipi 
     771               jj=1+floor(REAL(map%ptr(ib)-1)/REAL(ilendta)) 
     772               ji=map%ptr(ib)-(jj-1)*ilendta 
     773               DO ik = 1, ipk 
     774                  dta(ib,1,ik) =  dta_read(ji,jj,ik) 
     775               END DO 
     776            END DO 
     777         ENDIF 
     778 
    764779      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
    765780      ! Do we include something here to adjust barotropic velocities ! 
    766781      ! in case of a depth difference between bdy files and          ! 
    767782      ! bathymetry in the case ln_full_vel = .false. and jpk_bdy>0?  ! 
    768       ! [as the enveloping and parital cells could change H          ! 
     783      ! [as the enveloping and parital cells could change H]         ! 
    769784      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
     785 
     786      CASE DEFAULT   ;    
     787   
    770788 
    771789      IF( PRESENT(jpk_bdy) .AND. jpk_bdy>0 ) THEN       ! boundary data not on model grid: vertical interpolation 
     
    773791         SELECT CASE( igrd )                   
    774792            CASE(1) 
    775                WRITE(*,*) 'reading gdept'  
    776793               CALL iom_get ( num, jpdom_unknown, 'gdept', dta_read_z(1:ilendta,1:ipj,1:jpk_bdy) ) 
    777                WRITE(*,*) 'reading e3t'  
    778794               CALL iom_get ( num, jpdom_unknown, 'e3t', dta_read_dz(1:ilendta,1:ipj,1:jpk_bdy) ) 
    779                WRITE(*,*) 'finished reading'  
    780795            CASE(2)   
    781796               CALL iom_get ( num, jpdom_unknown, 'gdepu', dta_read_z(1:ilendta,1:ipj,1:jpk_bdy) ) 
     
    788803 
    789804#if defined key_bdy 
    790       WRITE(*,*) 'going in bdy interp'  
    791805         CALL fld_bdy_interp(dta_read, dta_read_z, dta_read_dz, map, jpk_bdy, igrd, ibdy, fv, dta, fvl) 
    792       WRITE(*,*) 'coming out bdy interp'  
    793806#endif 
    794807      ELSE ! boundary data assumed to be on model grid 
     808          
    795809         CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1:ipk), nrec )                     
    796810         IF (ipj==1) THEN ! we assume that this is an un-structured open boundary file 
Note: See TracChangeset for help on using the changeset viewer.