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_r11265_ASINTER-01_Guillaume_ABL1D/src/ABL – NEMO

source: NEMO/branches/2019/dev_r11265_ASINTER-01_Guillaume_ABL1D/src/ABL/sbcabl.F90 @ 11334

Last change on this file since 11334 was 11334, checked in by flemarie, 5 years ago

Follow-on to the ABL implementation (see ticket #2131)

  • Modification of diawri.F90 to output ABL variables with IOIPSL (i.e. without key_iomput)
  • Modification of ice_var_agg in icevar.F90 to allow the use of tm_su in ABL and bulk
  • Error handling in case ln_abl = .true. and ABL sources were not compiled
  • ABL now works with SI3 considering an average over ice categories
  • Update reference namelist to avoid runtime warnings due to nam_tide
  • Sanity checks with ORCA2 for the following configurations :

1 - ABL src + ln_blk = .true.
2 - ABL src + ln_abl = .true.
3 - no ABL src + ln_blk = .true.

All configurations are MPP reproducible and configurations 1 and 3 results are identical

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