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 5901 for branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 – NEMO

Ignore:
Timestamp:
2015-11-20T09:39:06+01:00 (8 years ago)
Author:
jamesharle
Message:

merging branch with head of the trunk

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90

    r5620 r5901  
    2929   USE daymod          ! calendar 
    3030   USE eosbn2          ! eq. of state, Brunt Vaisala frequency (eos     routine) 
    31    USE ldftra_oce      ! ocean active tracers: lateral physics 
     31   USE ldftra          ! lateral physics: ocean active tracers 
    3232   USE zdf_oce         ! ocean vertical physics 
    3333   USE phycst          ! physical constants 
    3434   USE dtatsd          ! data temperature and salinity   (dta_tsd routine) 
    3535   USE dtauvd          ! data: U & V current             (dta_uvd routine) 
    36    USE zpshde          ! partial step: hor. derivative (zps_hde routine) 
    37    USE eosbn2          ! equation of state            (eos bn2 routine) 
    3836   USE domvvl          ! varying vertical mesh 
    3937   USE dynspg_oce      ! pressure gradient schemes 
     
    7674      ! 
    7775 
    78       IF(lwp) WRITE(numout,*) ' ' 
     76      IF(lwp) WRITE(numout,*) 
    7977      IF(lwp) WRITE(numout,*) 'istate_ini : Initialization of the dynamics and tracers' 
    8078      IF(lwp) WRITE(numout,*) '~~~~~~~~~~' 
    8179 
    82       CALL dta_tsd_init                       ! Initialisation of T & S input data 
    83       IF( lk_c1d ) CALL dta_uvd_init          ! Initialization of U & V input data 
     80                     CALL dta_tsd_init        ! Initialisation of T & S input data 
     81      IF( lk_c1d )   CALL dta_uvd_init        ! Initialization of U & V input data 
    8482 
    8583      rhd  (:,:,:  ) = 0._wp   ;   rhop (:,:,:  ) = 0._wp      ! set one for all to 0 at level jpk 
     
    103101         ub   (:,:,:) = 0._wp   ;   un   (:,:,:) = 0._wp 
    104102         vb   (:,:,:) = 0._wp   ;   vn   (:,:,:) = 0._wp   
    105          rotb (:,:,:) = 0._wp   ;   rotn (:,:,:) = 0._wp 
    106          hdivb(:,:,:) = 0._wp   ;   hdivn(:,:,:) = 0._wp 
     103                                    hdivn(:,:,:) = 0._wp 
    107104         ! 
    108105         IF( cp_cfg == 'eel' ) THEN 
     
    119116            ENDIF 
    120117            IF ( ln_uvd_init .AND. lk_c1d ) THEN ! read 3D U and V data at nit000 
    121                CALL wrk_alloc( jpi, jpj, jpk, 2, zuvd ) 
     118               CALL wrk_alloc( jpi,jpj,jpk,2,  zuvd ) 
    122119               CALL dta_uvd( nit000, zuvd ) 
    123120               ub(:,:,:) = zuvd(:,:,:,1) ;  un(:,:,:) = ub(:,:,:) 
    124121               vb(:,:,:) = zuvd(:,:,:,2) ;  vn(:,:,:) = vb(:,:,:) 
    125                CALL wrk_dealloc( jpi, jpj, jpk, 2, zuvd ) 
     122               CALL wrk_dealloc( jpi,jpj,jpk,2,  zuvd ) 
    126123            ENDIF 
    127124         ENDIF 
    128          ! 
    129          CALL eos( tsb, rhd, rhop, gdept_0(:,:,:) )        ! before potential and in situ densities 
    130 #if ! defined key_c1d 
    131          IF( ln_zps .AND. .NOT. ln_isfcav)                                 & 
    132             &            CALL zps_hde    ( nit000, jpts, tsb, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient 
    133             &                                            rhd, gru , grv    )  ! of t, s, rd at the last ocean level 
    134          IF( ln_zps .AND.       ln_isfcav)                                 & 
    135             &            CALL zps_hde_isf( nit000, jpts, tsb, gtsu, gtsv,  &    ! Partial steps for top cell (ISF) 
    136             &                                            rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   & 
    137             &                                     gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level 
    138 #endif 
    139125         !    
    140126         ! - ML - sshn could be modified by istate_eel, so that initialization of fse3t_b is done here 
     
    142128            DO jk = 1, jpk 
    143129               fse3t_b(:,:,jk) = fse3t_n(:,:,jk) 
    144             ENDDO 
     130            END DO 
    145131         ENDIF 
    146132         !  
     
    163149      ! 
    164150      ! 
    165       un_b(:,:) = 0._wp ; vn_b(:,:) = 0._wp 
    166       ub_b(:,:) = 0._wp ; vb_b(:,:) = 0._wp 
     151      un_b(:,:) = 0._wp   ;  vn_b(:,:) = 0._wp 
     152      ub_b(:,:) = 0._wp   ;  vb_b(:,:) = 0._wp 
    167153      ! 
    168154      DO jk = 1, jpkm1 
     
    202188      !! References :  Philander ??? 
    203189      !!---------------------------------------------------------------------- 
    204       INTEGER  :: ji, jj, jk 
    205       REAL(wp) ::   zsal = 35.50 
     190      INTEGER  ::   ji, jj, jk 
     191      REAL(wp) ::   zsal = 35.50_wp 
    206192      !!---------------------------------------------------------------------- 
    207193      ! 
     
    210196      IF(lwp) WRITE(numout,*) '~~~~~~~~~~   and constant salinity (',zsal,' psu)' 
    211197      ! 
     198#if ! defined key_istate_fixed 
    212199      DO jk = 1, jpk 
    213200         tsn(:,:,jk,jp_tem) = (  ( ( 7.5 - 0. * ABS( gphit(:,:) )/30. ) * ( 1.-TANH((fsdept(:,:,jk)-80.)/30.) )   & 
     
    215202         tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem) 
    216203      END DO 
     204#else 
     205      tsn(:,:,:,jp_tem) = ztem * tmask(:,:,:) 
     206      tsb(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) 
     207#endif 
    217208      tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:) 
    218209      tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) 
     210 
    219211      ! 
    220212   END SUBROUTINE istate_t_s 
     
    233225      !!                and relative vorticity fields 
    234226      !!---------------------------------------------------------------------- 
    235       USE divcur     ! hor. divergence & rel. vorticity      (div_cur routine) 
     227      USE divhor     ! hor. divergence      (div_hor routine) 
    236228      USE iom 
    237229      ! 
     
    282274            tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) 
    283275            ! 
    284             ! set the dynamics: U,V, hdiv, rot (and ssh if necessary) 
     276            ! set the dynamics: U,V, hdiv (and ssh if necessary) 
    285277            ! ---------------- 
    286278            ! Start EEL5 configuration with barotropic geostrophic velocities  
     
    318310            ENDIF 
    319311            ! 
    320             CALL div_cur( nit000 )                  ! horizontal divergence and relative vorticity (curl) 
     312!!gm  Check  here call to div_hor should not be necessary 
     313!!gm         div_hor call runoffs  not sure they are defined at that level 
     314            CALL div_hor( nit000 )                  ! horizontal divergence and relative vorticity (curl) 
    321315            ! N.B. the vertical velocity will be computed from the horizontal divergence field 
    322316            ! in istate by a call to wzv routine 
     
    371365      !! 
    372366      !! ** Method  : - set temprature field 
    373       !!              - set salinity field 
     367      !!              - set salinity   field 
    374368      !!---------------------------------------------------------------------- 
    375369      INTEGER :: ji, jj, jk  ! dummy loop indices 
     
    458452      !!---------------------------------------------------------------------- 
    459453      USE dynspg          ! surface pressure gradient             (dyn_spg routine) 
    460       USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine) 
     454      USE divhor          ! hor. divergence                       (div_hor routine) 
    461455      USE lbclnk          ! ocean lateral boundary condition (or mpp link) 
    462456      ! 
     
    467461      !!---------------------------------------------------------------------- 
    468462      ! 
    469       CALL wrk_alloc( jpi, jpj, jpk, zprn) 
     463      CALL wrk_alloc( jpi,jpj,jpk,  zprn) 
    470464      ! 
    471465      IF(lwp) WRITE(numout,*)  
     
    544538      indic = 0 
    545539      CALL dyn_spg( nit000, indic )       ! surface pressure gradient 
    546  
     540      ! 
    547541      ! the new velocity is ua*rdt 
    548  
     542      ! 
    549543      CALL lbc_lnk( ua, 'U', -1. ) 
    550544      CALL lbc_lnk( va, 'V', -1. ) 
     
    556550      un(:,:,:) = ub(:,:,:) 
    557551      vn(:,:,:) = vb(:,:,:) 
    558         
    559       ! Compute the divergence and curl 
    560  
    561       CALL div_cur( nit000 )            ! now horizontal divergence and curl 
    562  
    563       hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value 
    564       rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value 
    565       ! 
    566       CALL wrk_dealloc( jpi, jpj, jpk, zprn) 
     552      ! 
     553!!gm  Check  here call to div_hor should not be necessary 
     554!!gm         div_hor call runoffs  not sure they are defined at that level 
     555      CALL div_hor( nit000 )            ! now horizontal divergence 
     556      ! 
     557      CALL wrk_dealloc( jpi,jpj,jpk,   zprn) 
    567558      ! 
    568559   END SUBROUTINE istate_uvg 
Note: See TracChangeset for help on using the changeset viewer.