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 11001 for NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/step.F90 – NEMO

Ignore:
Timestamp:
2019-05-20T14:19:17+02:00 (5 years ago)
Author:
davestorkey
Message:

2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps : C1D and step.F90. Passes SETTE. Compiles with key_c1d.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/step.F90

    r10998 r11001  
    125125      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    126126      IF( ln_sto_eos ) CALL sto_par( kstp )          ! Stochastic parameters 
    127       IF( ln_sto_eos ) CALL sto_pts( tsn  )          ! Random T/S fluctuations 
     127      IF( ln_sto_eos ) CALL sto_pts( ts(:,:,:,:,Nnn)  )          ! Random T/S fluctuations 
    128128 
    129129      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    131131      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    132132      !  THERMODYNAMICS 
    133                          CALL eos_rab( tsb, rab_b, Nnn )       ! before local thermal/haline expension ratio at T-points 
    134                          CALL eos_rab( tsn, rab_n, Nnn )       ! now    local thermal/haline expension ratio at T-points 
    135                          CALL bn2    ( tsb, rab_b, rn2b, Nnn ) ! before Brunt-Vaisala frequency 
    136                          CALL bn2    ( tsn, rab_n, rn2, Nnn  ) ! now    Brunt-Vaisala frequency 
     133                         CALL eos_rab( ts(:,:,:,:,Nbb), rab_b, Nnn )       ! before local thermal/haline expension ratio at T-points 
     134                         CALL eos_rab( ts(:,:,:,:,Nnn), rab_n, Nnn )       ! now    local thermal/haline expension ratio at T-points 
     135                         CALL bn2    ( ts(:,:,:,:,Nbb), rab_b, rn2b, Nnn ) ! before Brunt-Vaisala frequency 
     136                         CALL bn2    ( ts(:,:,:,:,Nnn), rab_n, rn2, Nnn  ) ! now    Brunt-Vaisala frequency 
    137137 
    138138      !  VERTICAL PHYSICS 
     
    142142      ! 
    143143      IF( l_ldfslp ) THEN                             ! slope of lateral mixing 
    144                          CALL eos( tsb, rhd, gdept_0(:,:,:) )               ! before in situ density 
     144                         CALL eos( ts(:,:,:,:,Nbb), rhd, gdept_0(:,:,:) )               ! before in situ density 
    145145 
    146146         IF( ln_zps .AND. .NOT. ln_isfcav)                                    & 
    147             &            CALL zps_hde    ( kstp, Nnn, jpts, tsb, gtsu, gtsv,  &  ! Partial steps: before horizontal gradient 
     147            &            CALL zps_hde    ( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv,  &  ! Partial steps: before horizontal gradient 
    148148            &                                          rhd, gru , grv    )       ! of t, s, rd at the last ocean level 
    149149 
    150150         IF( ln_zps .AND.       ln_isfcav)                                                & 
    151             &            CALL zps_hde_isf( kstp, Nnn, jpts, tsb, gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF) 
     151            &            CALL zps_hde_isf( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF) 
    152152            &                                          rhd, gru , grv , grui, grvi   )       ! of t, s, rd at the first ocean level 
    153153         IF( ln_traldf_triad ) THEN  
     
    169169                            CALL wzv           ( kstp, Nbb, Nnn, ww,  Naa )    ! now cross-level velocity  
    170170      IF( ln_zad_Aimp )     CALL wAimp         ( kstp,      Nnn           )  ! Adaptive-implicit vertical advection partitioning 
    171                             CALL eos    ( tsn, rhd, rhop, gdept_n(:,:,:) )  ! now in situ density for hpg computation 
     171                            CALL eos    ( ts(:,:,:,:,Nnn), rhd, rhop, gdept(:,:,:,Nnn) )  ! now in situ density for hpg computation 
    172172                             
    173173!!jc: fs simplification 
     
    176176!!                                         with previous versions using split-explicit free surface           
    177177            IF( ln_zps .AND. .NOT. ln_isfcav )                                    & 
    178                &            CALL zps_hde    ( kstp, Nnn, jpts, tsn, gtsu, gtsv,   &  ! Partial steps: before horizontal gradient 
     178               &            CALL zps_hde    ( kstp, Nnn, jpts, ts(:,:,:,:,Nnn), gtsu, gtsv,   &  ! Partial steps: before horizontal gradient 
    179179               &                                          rhd, gru , grv     )       ! of t, s, rd at the last ocean level 
    180180            IF( ln_zps .AND.       ln_isfcav )                                               & 
    181                &            CALL zps_hde_isf( kstp, Nnn, jpts, tsn, gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF) 
     181               &            CALL zps_hde_isf( kstp, Nnn, jpts, ts(:,:,:,:,Nnn), gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF) 
    182182               &                                          rhd, gru , grv , grui, grvi   )       ! of t, s, rd at the first ocean level 
    183183!!jc: fs simplification 
    184184                             
    185                          ua(:,:,:) = 0._wp            ! set dynamics trends to zero 
    186                          va(:,:,:) = 0._wp 
     185                         uu(:,:,:,Nrhs) = 0._wp            ! set dynamics trends to zero 
     186                         vv(:,:,:,Nrhs) = 0._wp 
    187187 
    188188      IF(  lk_asminc .AND. ln_asmiau .AND. ln_dyninc )   & 
     
    200200                         CALL dyn_spg( kstp, Nbb, Nnn, Nrhs, uu, vv, ssh, uu_b, vv_b, Naa )  ! surface pressure gradient 
    201201 
    202                                                       ! With split-explicit free surface, since now transports have been updated and ssha as well 
     202                                                      ! With split-explicit free surface, since now transports have been updated and ssh(:,:,Nrhs) as well 
    203203      IF( ln_dynspg_ts ) THEN                         ! vertical scale factors and vertical velocity need to be updated 
    204204                            CALL div_hor       ( kstp, Nbb, Nnn )    ! Horizontal divergence  (2nd call in time-split case) 
     
    238238      ! Active tracers                               
    239239      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    240                          tsa(:,:,:,:) = 0._wp         ! set tracer trends to zero 
     240                         ts(:,:,:,:,Nrhs) = 0._wp         ! set tracer trends to zero 
    241241 
    242242      IF(  lk_asminc .AND. ln_asmiau .AND. & 
    243          & ln_trainc )   CALL tra_asm_inc   ( kstp, Nbb, Nnn, ts, Nrhs )  ! apply tracer assimilation increment 
    244                          CALL tra_sbc       ( kstp,      Nnn, ts, Nrhs )  ! surface boundary condition 
    245       IF( ln_traqsr  )   CALL tra_qsr       ( kstp,      Nnn, ts, Nrhs )  ! penetrative solar radiation qsr 
    246       IF( ln_trabbc  )   CALL tra_bbc       ( kstp,      Nnn, ts, Nrhs )  ! bottom heat flux 
    247       IF( ln_trabbl  )   CALL tra_bbl       ( kstp, Nbb, Nnn, ts, Nrhs )  ! advective (and/or diffusive) bottom boundary layer scheme 
    248       IF( ln_tradmp  )   CALL tra_dmp       ( kstp, Nbb, Nnn, ts, Nrhs )  ! internal damping trends 
    249       IF( ln_bdy     )   CALL bdy_tra_dmp   ( kstp, Nbb,      ts, Nrhs )  ! bdy damping trends 
     243         & ln_trainc )   CALL tra_asm_inc( kstp, Nbb, Nnn, ts, Nrhs )  ! apply tracer assimilation increment 
     244                         CALL tra_sbc    ( kstp,      Nnn, ts, Nrhs )  ! surface boundary condition 
     245      IF( ln_traqsr  )   CALL tra_qsr    ( kstp,      Nnn, ts, Nrhs )  ! penetrative solar radiation qsr 
     246      IF( ln_trabbc  )   CALL tra_bbc    ( kstp,      Nnn, ts, Nrhs )  ! bottom heat flux 
     247      IF( ln_trabbl  )   CALL tra_bbl    ( kstp, Nbb, Nnn, ts, Nrhs )  ! advective (and/or diffusive) bottom boundary layer scheme 
     248      IF( ln_tradmp  )   CALL tra_dmp    ( kstp, Nbb, Nnn, ts, Nrhs )  ! internal damping trends 
     249      IF( ln_bdy     )   CALL bdy_tra_dmp( kstp, Nbb,      ts, Nrhs )  ! bdy damping trends 
    250250#if defined key_agrif 
    251251      IF(.NOT. Agrif_Root())  &  
    252252               &         CALL Agrif_Sponge_tra        ! tracers sponge 
    253253#endif 
    254                          CALL tra_adv( kstp, Nbb, Nnn      , ts, Nrhs )  ! hor. + vert. advection  ==> RHS 
    255       IF( ln_zdfosm  )   CALL tra_osm( kstp,      Nnn      , ts, Nrhs )  ! OSMOSIS non-local tracer fluxes ==> RHS 
     254                         CALL tra_adv    ( kstp, Nbb, Nnn, ts, Nrhs )  ! hor. + vert. advection ==> RHS 
     255      IF( ln_zdfosm  )   CALL tra_osm    ( kstp, Nnn     , ts, Nrhs )  ! OSMOSIS non-local tracer fluxes ==> RHS 
    256256      IF( lrst_oce .AND. ln_zdfosm ) & 
    257            &             CALL osm_rst( kstp,      Nnn, 'WRITE' )         ! write OSMOSIS outputs + wn (so must do here) to restarts 
    258                          CALL tra_ldf( kstp, Nbb, Nnn      , ts, Nrhs )  ! lateral mixing 
     257           &             CALL osm_rst    ( kstp, Nnn, 'WRITE' )        ! write OSMOSIS outputs + ww (so must do here) to restarts 
     258                         CALL tra_ldf    ( kstp, Nbb, Nnn, ts, Nrhs )  ! lateral mixing 
    259259 
    260260!!gm : why CALL to dia_ptr has been moved here??? (use trends info?) 
    261261      IF( ln_diaptr  )   CALL dia_ptr( Nnn )                 ! Poleward adv/ldf TRansports diagnostics 
    262262!!gm 
    263                          CALL tra_zdf( kstp, Nbb, Nnn, Nrhs, ts, Naa  )  ! vert. mixing & after tracer   ==> after 
    264       IF( ln_zdfnpc  )   CALL tra_npc( kstp,      Nnn, Nrhs, ts, Naa  )  ! update after fields by non-penetrative convection 
     263                         CALL tra_zdf    ( kstp, Nbb, Nnn, Nrhs, ts, Naa  )  ! vert. mixing & after tracer  ==> after 
     264      IF( ln_zdfnpc  )   CALL tra_npc    ( kstp,      Nnn, Nrhs, ts, Naa  )  ! update after fields by non-penetrative convection 
    265265 
    266266      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    280280!!    place. 
    281281!!  
    282 !!jc2: dynnxt must be the latest call. e3t_b are indeed updated in that routine 
     282!!jc2: dynnxt must be the latest call. e3t(:,:,:,Nbb) are indeed updated in that routine 
    283283                         CALL tra_nxt       ( kstp, Nbb, Nnn, Nrhs, Naa )  ! finalize (bcs) tracer fields at next time step and swap 
    284284                         CALL dyn_nxt       ( kstp, Nbb, Nnn, Naa  )  ! finalize (bcs) velocities at next time step and swap (always called after tra_nxt) 
Note: See TracChangeset for help on using the changeset viewer.