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.
sbcabl.F90 in NEMO/branches/2019/dev_r11943_MERGE_2019/src/ABL – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/ABL/sbcabl.F90 @ 12199

Last change on this file since 12199 was 12199, checked in by laurent, 4 years ago

Catches up with SBCBLK bugfixes of branch "dev_r12072_MERGE_OPTION2_2019"

File size: 20.2 KB
RevLine 
[11305]1MODULE sbcabl
2   !!======================================================================
3   !!                       ***  MODULE  sbcabl  ***
[12199]4   !! Ocean forcing:  momentum, heat and freshwater flux formulation
5   !!                 derived from an ABL model
[11305]6   !!=====================================================================
7   !! History :  4.0  !  2019-03  (F. Lemarié & G. Samson)  Original code
8   !!----------------------------------------------------------------------
9
10   !!----------------------------------------------------------------------
[12199]11   !!   sbc_abl_init  : Initialization of ABL model based on namelist options
12   !!   sbc_abl       : driver for the computation of momentum, heat and freshwater
[11305]13   !!                   fluxes over ocean via the ABL model
14   !!----------------------------------------------------------------------
15   USE abl            ! ABL
16   USE par_abl        ! abl parameters
17   USE ablmod
[11858]18   USE ablrst
[11305]19
20   USE phycst         ! physical constants
21   USE fldread        ! read input fields
22   USE sbc_oce        ! Surface boundary condition: ocean fields
23   USE sbcblk         ! Surface boundary condition: bulk formulae
[12021]24   USE sbcblk_phy     ! Surface boundary condition: bulk formulae
[11305]25   USE dom_oce, ONLY  : tmask
26   !
27   USE iom            ! I/O manager library
28   USE in_out_manager ! I/O manager
29   USE lib_mpp        ! distribued memory computing library
30   USE lib_fortran    ! to use key_nosignedzero
31   USE timing         ! Timing
32   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
33   USE prtctl         ! Print control
[11334]34#if defined key_si3
35   USE ice    , ONLY : u_ice, v_ice, tm_su, ato_i      ! ato_i = total open water fractional area
36   USE sbc_ice, ONLY : wndm_ice, utau_ice, vtau_ice
[12199]37#endif
[11334]38#if ! defined key_iomput
39   USE diawri    , ONLY : dia_wri_alloc_abl
40#endif
[11305]41   IMPLICIT NONE
42   PRIVATE
43
44   PUBLIC   sbc_abl_init       ! routine called in sbcmod module
45   PUBLIC   sbc_abl            ! routine called in sbcmod module
46
47   !! * Substitutions
48#  include "vectopt_loop_substitute.h90"
49   !!----------------------------------------------------------------------
50   !! NEMO/OPA 3.7 , NEMO-consortium (2014)
51   !! $Id: sbcabl.F90 6416 2016-04-01 12:22:17Z clem $
52   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
53   !!----------------------------------------------------------------------
54CONTAINS
55
56   SUBROUTINE sbc_abl_init
57      !!---------------------------------------------------------------------
58      !!                    ***  ROUTINE sbc_abl_init  ***
59      !!
60      !! ** Purposes :   - read namelist section namsbc_abl
[12199]61      !!                 - initialize and check parameter values
62      !!                 - initialize variables of ABL model
[11305]63      !!
64      !!----------------------------------------------------------------------
65      INTEGER            ::   ji, jj, jk, jbak, jbak_dta       ! dummy loop indices
66      INTEGER            ::   ios, ierror, ioptio              ! Local integer
67      INTEGER            ::   inum, indims, idimsz(4), id
68      CHARACTER(len=100) ::   cn_dir, cn_dom                   ! Atmospheric grid directory
69      REAL(wp)           ::   zcff,zcff1
70      LOGICAL            ::   lluldl
[11858]71      NAMELIST/namsbc_abl/ cn_dir, cn_dom, cn_ablrst_in, cn_ablrst_out,           &
72         &                 cn_ablrst_indir, cn_ablrst_outdir,                     &
73         &                 ln_hpgls_frc, ln_geos_winds, nn_dyn_restore,           &
[11305]74         &                 rn_ldyn_min , rn_ldyn_max, rn_ltra_min, rn_ltra_max,   &
[11322]75         &                 nn_amxl, rn_cm, rn_ct, rn_ce, rn_ceps, rn_Rod, rn_Ric, &
[12199]76         &                 ln_smth_pblh
[11305]77      !!---------------------------------------------------------------------
78
[12199]79      REWIND( numnam_ref )              ! Namelist namsbc_abl in reference namelist : ABL parameters
[11305]80      READ  ( numnam_ref, namsbc_abl, IOSTAT = ios, ERR = 901 )
[11348]81901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_abl in reference namelist' )
[11305]82      !
[12199]83      REWIND( numnam_cfg )              ! Namelist namsbc_abl in configuration namelist : ABL parameters
[11305]84      READ  ( numnam_cfg, namsbc_abl, IOSTAT = ios, ERR = 902 )
[11348]85902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_abl in configuration namelist' )
[11305]86      !
87      IF(lwm) WRITE( numond, namsbc_abl )
88      !
[12199]89      ! Check ABL mixing length option
[11305]90      IF( nn_amxl  < 0   .OR.  nn_amxl  > 2 )   &
[12199]91         &                 CALL ctl_stop( 'abl_init : bad flag, nn_amxl must be  0, 1 or 2 ' )
[11305]92      !
[12199]93      ! Check ABL dyn restore option
[11305]94      IF( nn_dyn_restore  < 0   .OR.  nn_dyn_restore  > 2 )   &
[12199]95         &                 CALL ctl_stop( 'abl_init : bad flag, nn_dyn_restore must be  0, 1 or 2 ' )
[12063]96
[11305]97      !!---------------------------------------------------------------------
98      !! Control prints
[12199]99      !!---------------------------------------------------------------------
[11305]100      IF(lwp) THEN                              ! Control print (other namelist variable)
101         WRITE(numout,*)
102         WRITE(numout,*) '      ABL -- cn_dir         = ', cn_dir
103         WRITE(numout,*) '      ABL -- cn_dom         = ', cn_dom
104         IF( ln_hpgls_frc ) THEN
105            WRITE(numout,*) '      ABL -- winds forced by large-scale pressure gradient'
106            IF(ln_geos_winds) THEN
107               ln_geos_winds = .FALSE.
[12199]108               WRITE(numout,*) '      ABL -- geostrophic guide disabled (not compatible with ln_hpgls_frc = .T.)'
[11305]109            END IF
110         ELSE IF( ln_geos_winds ) THEN
[12199]111            WRITE(numout,*) '      ABL -- winds forced by geostrophic winds'
[11305]112         ELSE
[12199]113            WRITE(numout,*) '      ABL -- Geostrophic winds and large-scale pressure gradient are ignored'
[11305]114         END IF
115         !
116         SELECT CASE ( nn_dyn_restore )
[12199]117         CASE ( 0 )
118            WRITE(numout,*) '      ABL -- No restoring for ABL winds'
[11305]119         CASE ( 1 )
[12199]120            WRITE(numout,*) '      ABL -- Restoring of ABL winds only in the equatorial region '
[11305]121         CASE ( 2 )
[12199]122            WRITE(numout,*) '      ABL -- Restoring of ABL winds activated everywhere '
[11305]123         END SELECT
124         !
[12199]125         IF( ln_smth_pblh )  WRITE(numout,*) '      ABL -- Smoothing of PBL height is activated'
126         !
[11305]127      ENDIF
128
129      !!---------------------------------------------------------------------
130      !! Convert nudging coefficient from hours to 1/sec
[12199]131      !!---------------------------------------------------------------------
[11305]132      zcff        = 1._wp / 3600._wp
133      rn_ldyn_min = zcff / rn_ldyn_min
[12199]134      rn_ldyn_max = zcff / rn_ldyn_max
[11305]135      rn_ltra_min = zcff / rn_ltra_min
[12199]136      rn_ltra_max = zcff / rn_ltra_max
[11305]137
138      !!---------------------------------------------------------------------
[12199]139      !! ABL grid initialization
140      !!---------------------------------------------------------------------
141      CALL iom_open( TRIM(cn_dir)//TRIM(cn_dom), inum )
[11305]142      id     = iom_varid( inum, 'e3t_abl', kdimsz=idimsz, kndims=indims, lduld=lluldl )
[12199]143      jpka   = idimsz(indims - COUNT( (/lluldl/) ) )
[11305]144      jpkam1 = jpka - 1
145
146      IF( abl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'abl_init : unable to allocate arrays' )
147      CALL iom_get( inum, jpdom_unknown, 'e3t_abl', e3t_abl(:) )
148      CALL iom_get( inum, jpdom_unknown, 'e3w_abl', e3w_abl(:) )
149      CALL iom_get( inum, jpdom_unknown, 'ght_abl', ght_abl(:) )
150      CALL iom_get( inum, jpdom_unknown, 'ghw_abl', ghw_abl(:) )
151      CALL iom_close( inum )
[11348]152
[11334]153#if ! defined key_iomput
[12199]154      IF( dia_wri_alloc_abl()  /= 0 ) CALL ctl_stop( 'STOP', 'abl_init : unable to allocate arrays' )
[11334]155#endif
[11305]156
157      IF(lwp) THEN
158         WRITE(numout,*)
159         WRITE(numout,*) '    sbc_abl_init   : ABL Reference vertical grid'
[12199]160         WRITE(numout,*) '    ~~~~~~~'
[11305]161         WRITE(numout, "(9x,'  level ght_abl   ghw_abl   e3t_abl   e3w_abl  ')" )
162         WRITE(numout, "(10x, i4, 4f9.2)" ) ( jk, ght_abl(jk), ghw_abl(jk), e3t_abl(jk), e3w_abl(jk), jk = 1, jpka )
163      END IF
164
165      !!---------------------------------------------------------------------
[12199]166      !! Check TKE closure parameters
167      !!---------------------------------------------------------------------
[11305]168      rn_Sch  = rn_ce / rn_cm
169      mxl_min = (avm_bak / rn_cm) / sqrt( tke_min )
170
171      IF(lwp) THEN
172         WRITE(numout,*)
173         WRITE(numout,*) '    abl_zdf_tke   : ABL TKE turbulent closure'
[12199]174         WRITE(numout,*) '    ~~~~~~~~~~~'
[11305]175         IF(nn_amxl==0) WRITE(numout,*) 'Deardorff 80 length-scale '
176         IF(nn_amxl==1) WRITE(numout,*) 'length-scale based on the distance to the PBL height '
177         WRITE(numout,*) ' Minimum value of atmospheric TKE           = ',tke_min,' m^2 s^-2'
178         WRITE(numout,*) ' Minimum value of atmospheric mixing length = ',mxl_min,' m'
[12199]179         WRITE(numout,*) ' Constant for turbulent viscosity           = ',rn_Cm
180         WRITE(numout,*) ' Constant for turbulent diffusivity         = ',rn_Ct
181         WRITE(numout,*) ' Constant for Schmidt number                = ',rn_Sch
182         WRITE(numout,*) ' Constant for TKE dissipation               = ',rn_Ceps
[11305]183      END IF
184
185      !!-------------------------------------------------------------------------------------------
186      !! Compute parameters to build the vertical profile for the nudging term (used in abl_stp())
[12199]187      !!-------------------------------------------------------------------------------------------
[11305]188      zcff1       = 1._wp / ( jp_bmax - jp_bmin )**3
189      ! for active tracers
190      jp_alp3_tra = -2._wp * zcff1                       * ( rn_ltra_max - rn_ltra_min )
191      jp_alp2_tra =  3._wp * zcff1 * (jp_bmax + jp_bmin) * ( rn_ltra_max - rn_ltra_min )
192      jp_alp1_tra = -6._wp * zcff1 *  jp_bmax * jp_bmin  * ( rn_ltra_max - rn_ltra_min )
193      jp_alp0_tra = zcff1 * (  rn_ltra_max * jp_bmin*jp_bmin * (3._wp*jp_bmax - jp_bmin)      &
[12199]194         &                   - rn_ltra_min * jp_bmax*jp_bmax * (3._wp*jp_bmin - jp_bmax) )
[11305]195      ! for dynamics
196      jp_alp3_dyn = -2._wp * zcff1                       * ( rn_ldyn_max - rn_ldyn_min )
197      jp_alp2_dyn =  3._wp * zcff1 * (jp_bmax + jp_bmin) * ( rn_ldyn_max - rn_ldyn_min )
198      jp_alp1_dyn = -6._wp * zcff1 *  jp_bmax * jp_bmin  * ( rn_ldyn_max - rn_ldyn_min )
199      jp_alp0_dyn = zcff1 * (  rn_ldyn_max * jp_bmin*jp_bmin * (3._wp*jp_bmax - jp_bmin)      &
[12199]200         &                   - rn_ldyn_min * jp_bmax*jp_bmax * (3._wp*jp_bmin - jp_bmax) )
[11305]201
202      jp_pblh_min = ghw_abl(     4) / jp_bmin  !<-- at least 3 grid points at the bottom have value rn_ltra_min
203      jp_pblh_max = ghw_abl(jpka-3) / jp_bmax  !<-- at least 3 grid points at the top    have value rn_ltra_max
204
[11363]205      ! ABL timestep
206      rdt_abl = nn_fsbc * rdt
207
[11305]208      ! Check parameters for dynamics
[11363]209      zcff  = ( jp_alp3_dyn * jp_bmin**3 + jp_alp2_dyn * jp_bmin**2   &
210         &    + jp_alp1_dyn * jp_bmin    + jp_alp0_dyn ) * rdt_abl
211      zcff1 = ( jp_alp3_dyn * jp_bmax**3 + jp_alp2_dyn * jp_bmax**2   &
212         &    + jp_alp1_dyn * jp_bmax    + jp_alp0_dyn ) * rdt_abl
[11305]213      IF(lwp) THEN
214         IF(nn_dyn_restore > 0) THEN
[11363]215            WRITE(numout,*) ' ABL Minimum value for dynamics restoring = ',zcff
216            WRITE(numout,*) ' ABL Maximum value for dynamics restoring = ',zcff1
[11305]217            ! Check that restoring coefficients are between 0 and 1
[11586]218            IF( zcff1 - nn_fsbc > 0.001_wp .OR. zcff1 < 0._wp )   &
[11305]219               &                   CALL ctl_stop( 'abl_init : wrong value for rn_ldyn_max' )
[11586]220            IF( zcff  - nn_fsbc > 0.001_wp .OR. zcff  < 0._wp )   &
[11305]221               &                   CALL ctl_stop( 'abl_init : wrong value for rn_ldyn_min' )
[11363]222            IF( zcff  > zcff1                    )   &
[11586]223               &                   CALL ctl_stop( 'abl_init : rn_ldyn_max must be smaller than rn_ldyn_min' )
[11305]224         END IF
225      END IF
226
227      ! Check parameters for active tracers
[11363]228      zcff  = ( jp_alp3_tra * jp_bmin**3 + jp_alp2_tra * jp_bmin**2   &
229         &    + jp_alp1_tra * jp_bmin    + jp_alp0_tra ) * rdt_abl
230      zcff1 = ( jp_alp3_tra * jp_bmax**3 + jp_alp2_tra * jp_bmax**2   &
[12199]231         &    + jp_alp1_tra * jp_bmax    + jp_alp0_tra ) * rdt_abl
[11305]232      IF(lwp) THEN
[11363]233         WRITE(numout,*) ' ABL Minimum value for tracers restoring = ',zcff
234         WRITE(numout,*) ' ABL Maximum value for tracers restoring = ',zcff1
[11305]235         ! Check that restoring coefficients are between 0 and 1
[11586]236         IF( zcff1 - nn_fsbc > 0.001_wp .OR. zcff1 < 0._wp )   &
[11305]237            &                   CALL ctl_stop( 'abl_init : wrong value for rn_ltra_max' )
[11586]238         IF( zcff  - nn_fsbc > 0.001_wp .OR. zcff  < 0._wp )   &
[11305]239            &                   CALL ctl_stop( 'abl_init : wrong value for rn_ltra_min' )
[11363]240         IF( zcff  > zcff1                    )   &
[11586]241            &                   CALL ctl_stop( 'abl_init : rn_ltra_max must be smaller than rn_ltra_min' )
[11305]242      END IF
243
244      !!-------------------------------------------------------------------------------------------
245      !! Initialize Coriolis frequency, equatorial restoring and land/sea mask
[12199]246      !!-------------------------------------------------------------------------------------------
[11305]247      fft_abl(:,:) = 2._wp * omega * SIN( rad * gphit(:,:) )
248
249      ! Equatorial restoring
250      IF( nn_dyn_restore == 1 ) THEN
[12199]251         zcff         = 2._wp * omega * SIN( rad * 90._wp )   !++ fmax
[11363]252         rest_eq(:,:) = SIN( 0.5_wp*rpi*( (fft_abl(:,:) - zcff) / zcff ) )**8
253         !!GS: alternative shape
254         !rest_eq(:,:) = SIN( 0.5_wp*rpi*(zcff - ABS(ff_t(:,:))) / (zcff - 3.e-5) )**8
255         !WHERE(ABS(ff_t(:,:)).LE.3.e-5) rest_eq(:,:) = 1._wp
[11305]256      ELSE
257         rest_eq(:,:) = 1._wp
258      END IF
259      ! T-mask
260      msk_abl(:,:) = tmask(:,:,1)
261
[12199]262      !!-------------------------------------------------------------------------------------------
[11305]263
264      ! initialize 2D bulk fields AND 3D abl data
265      CALL sbc_blk_init
266
267      ! Initialize the time index for now time (nt_n) and after time (nt_a)
268      nt_n = 1 + MOD( nit000  , 2)
269      nt_a = 1 + MOD( nit000+1, 2)
270
[11858]271      ! initialize ABL from data or restart
272      IF( ln_rstart ) THEN
273         CALL abl_rst_read
274      ELSE
275         CALL fld_read( nit000, nn_fsbc, sf ) ! input fields provided at the first time-step
[11305]276
[12199]277         u_abl(:,:,:,nt_n      ) = sf(jp_wndi)%fnow(:,:,:)
278         v_abl(:,:,:,nt_n      ) = sf(jp_wndj)%fnow(:,:,:)
[11858]279         tq_abl(:,:,:,nt_n,jp_ta) = sf(jp_tair)%fnow(:,:,:)
280         tq_abl(:,:,:,nt_n,jp_qa) = sf(jp_humi)%fnow(:,:,:)
[12199]281
[11858]282         tke_abl(:,:,:,nt_n     ) = tke_min
283         avm_abl(:,:,:          ) = avm_bak
284         avt_abl(:,:,:          ) = avt_bak
[12199]285         mxl_abl(:,:,:          ) = mxl_min
286         pblh   (:,:            ) = ghw_abl( 3 )  !<-- assume that the pbl contains 3 grid points
[11858]287         u_abl  (:,:,:,nt_a     ) = 0._wp
288         v_abl  (:,:,:,nt_a     ) = 0._wp
289         tq_abl (:,:,:,nt_a,:   ) = 0._wp
290         tke_abl(:,:,:,nt_a     ) = 0._wp
291      ENDIF
[12021]292
293      rhoa(:,:) = rho_air( tq_abl(:,:,2,nt_n,jp_ta), tq_abl(:,:,2,nt_n,jp_qa), sf(jp_slp)%fnow(:,:,1) ) !!GS: rhoa must be (re)computed here here to avoid division by zero in blk_ice_1 (TBI)
[12199]294
[11305]295   END SUBROUTINE sbc_abl_init
296
297
298   SUBROUTINE sbc_abl( kt )
299      !!---------------------------------------------------------------------
300      !!                     ***  ROUTINE sbc_abl  ***
301      !!
302      !! ** Purpose :   provide the momentum, heat and freshwater fluxes at
[12199]303      !!      the ocean surface from an ABL calculation at each oceanic time step
[11305]304      !!
[12199]305      !! ** Method  :
[11305]306      !!              - Pre-compute part of turbulent fluxes in blk_oce_1
[12199]307      !!              - Perform 1 time-step of the ABL model
[11305]308      !!              - Finalize flux computation in blk_oce_2
309      !!
310      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
311      !!              - vtau    : j-component of the stress at V-point  (N/m2)
312      !!              - taum    : Wind stress module at T-point         (N/m2)
313      !!              - wndm    : Wind speed module at T-point          (m/s)
314      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
315      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
316      !!              - emp     : evaporation minus precipitation       (kg/m2/s)
317      !!
318      !!---------------------------------------------------------------------
319      INTEGER ,         INTENT(in) ::   kt   ! ocean time step
320      !!
321      REAL(wp), DIMENSION(jpi,jpj) ::   zssq, zcd_du, zsen, zevp
[11334]322#if defined key_si3
323      REAL(wp), DIMENSION(jpi,jpj) ::   zssqi, zcd_dui, zseni, zevpi
[12199]324#endif
[11305]325      INTEGER                      ::   jbak, jbak_dta, ji, jj
326      !!---------------------------------------------------------------------
327      !
328      !!-------------------------------------------------------------------------------------------
329      !! 1 - Read Atmospheric 3D data for large-scale forcing
330      !!-------------------------------------------------------------------------------------------
[12199]331
[11305]332      CALL fld_read( kt, nn_fsbc, sf )             ! input fields provided at the current time-step
[12199]333
[11305]334      !!-------------------------------------------------------------------------------------------
335      !! 2 - Compute Cd x ||U||, Ch x ||U||, Ce x ||U||, and SSQ using now fields
[12199]336      !!-------------------------------------------------------------------------------------------
[11334]337
[11348]338      CALL blk_oce_1( kt,  u_abl(:,:,2,nt_n      ),  v_abl(:,:,2,nt_n      ),   &   !   <<= in
339         &                tq_abl(:,:,2,nt_n,jp_ta), tq_abl(:,:,2,nt_n,jp_qa),   &   !   <<= in
340         &                sf(jp_slp )%fnow(:,:,1) , sst_m, ssu_m, ssv_m     ,   &   !   <<= in
[12015]341         &                sf(jp_qsr )%fnow(:,:,1) , sf(jp_qlw )%fnow(:,:,1) ,   &   !   <<= in
[12199]342         &                tsk_m, zssq, zcd_du, zsen, zevp                       )   !   =>> out
343
[11334]344#if defined key_si3
[11360]345      CALL blk_ice_1(  u_abl(:,:,2,nt_n      ),  v_abl(:,:,2,nt_n      ),    &   !   <<= in
346         &            tq_abl(:,:,2,nt_n,jp_ta), tq_abl(:,:,2,nt_n,jp_qa),    &   !   <<= in
347         &            sf(jp_slp)%fnow(:,:,1)  ,  u_ice, v_ice, tm_su    ,    &   !   <<= in
348         &            pseni=zseni, pevpi=zevpi, pssqi=zssqi, pcd_dui=zcd_dui )   !   <<= out
[12199]349#endif
[11334]350
[11305]351      !!-------------------------------------------------------------------------------------------
352      !! 3 - Advance ABL variables from now (n) to after (n+1)
[12199]353      !!-------------------------------------------------------------------------------------------
354
355      CALL abl_stp( kt, tsk_m, ssu_m, ssv_m, zssq,                          &   !   <<= in
[11348]356         &              sf(jp_wndi)%fnow(:,:,:), sf(jp_wndj)%fnow(:,:,:),   &   !   <<= in
357         &              sf(jp_tair)%fnow(:,:,:), sf(jp_humi)%fnow(:,:,:),   &   !   <<= in
358         &              sf(jp_slp )%fnow(:,:,1),                            &   !   <<= in
359         &              sf(jp_hpgi)%fnow(:,:,:), sf(jp_hpgj)%fnow(:,:,:),   &   !   <<= in
360         &              zcd_du, zsen, zevp,                                 &   !   <=> in/out
361         &              wndm, utau, vtau, taum                              &   !   =>> out
[12199]362#if defined key_si3
[11348]363         &            , tm_su, u_ice, v_ice, zssqi, zcd_dui                 &   !   <<= in
364         &            , zseni, zevpi, wndm_ice, ato_i                       &   !   <<= in
365         &            , utau_ice, vtau_ice                                  &   !   =>> out
[12199]366#endif
[11348]367         &                                                                  )
[11305]368      !!-------------------------------------------------------------------------------------------
[12199]369      !! 4 - Finalize flux computation using ABL variables at (n+1), nt_n corresponds to (n+1) since
[11305]370      !!                                                                time swap is done in abl_stp
[12199]371      !!-------------------------------------------------------------------------------------------
[11305]372
[11360]373      CALL blk_oce_2( tq_abl(:,:,2,nt_n,jp_ta),                            &
374         &            sf(jp_qsr )%fnow(:,:,1) , sf(jp_qlw )%fnow(:,:,1),   &
375         &            sf(jp_prec)%fnow(:,:,1) , sf(jp_snow)%fnow(:,:,1),   &
[12199]376         &            tsk_m, zsen, zevp                                )
[11305]377
[12199]378      CALL abl_rst_opn( kt )                       ! Open abl restart file (if necessary)
379      IF( lrst_abl ) CALL abl_rst_write( kt )      ! -- abl restart file
[11858]380
[11334]381#if defined key_si3
[11348]382      ! Avoid a USE abl in icesbc module
[11334]383      sf(jp_tair)%fnow(:,:,1) = tq_abl(:,:,2,nt_n,jp_ta);  sf(jp_humi)%fnow(:,:,1) = tq_abl(:,:,2,nt_n,jp_qa)
[12199]384#endif
[11334]385
[11305]386   END SUBROUTINE sbc_abl
387
388   !!======================================================================
389END MODULE sbcabl
Note: See TracBrowser for help on using the repository browser.