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.
sbcblk_core.F90 in trunk/NEMO/OPA_SRC/SBC – NEMO

source: trunk/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 1037

Last change on this file since 1037 was 1037, checked in by ctlod, 16 years ago

trunk: replace freeze(:,:) variable with fr_i(:,:), use the tfreez function defined in eosbn2.F90 and remove the useless ocfzpt.F90 module, see ticket: #177

File size: 45.2 KB
Line 
1MODULE sbcblk_core
2   !!======================================================================
3   !!                       ***  MODULE  sbcblk_core  ***
4   !! Ocean forcing:  momentum, heat and freshwater flux formulation
5   !!=====================================================================
6   !! History :  1.0   !  04-08  (U. Schweckendiek)  Original code
7   !!            2.0   !  05-04  (L. Brodeau, A.M. Treguier) additions:
8   !!                            -  new bulk routine for efficiency
9   !!                            -  WINDS ARE NOW ASSUMED TO BE AT T POINTS in input files !!!!
10   !!                            -  file names and file characteristics in namelist
11   !!                            -  Implement reading of 6-hourly fields   
12   !!            3.0   !  06-06  (G. Madec) sbc rewritting
13   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   sbc_blk_core  : bulk formulation as ocean surface boundary condition
17   !!                   (forced mode, CORE bulk formulea)
18   !!   blk_oce_core  : ocean: computes momentum, heat and freshwater fluxes
19   !!   blk_ice_core  : ice  : computes momentum, heat and freshwater fluxes
20   !!   turb_core     : computes the CORE turbulent transfer coefficients
21   !!----------------------------------------------------------------------
22   USE oce             ! ocean dynamics and tracers
23   USE dom_oce         ! ocean space and time domain
24   USE phycst          ! physical constants
25   USE daymod          ! calendar
26   USE fldread         ! read input fields
27   USE sbc_oce         ! Surface boundary condition: ocean fields
28   USE iom             ! I/O manager library
29   USE in_out_manager  ! I/O manager
30   USE lib_mpp         ! distribued memory computing library
31   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
32   USE prtctl          ! Print control
33#if defined key_lim3
34   USE ice_oce         ! For ice surface temperature
35#endif
36
37
38   IMPLICIT NONE
39   PRIVATE
40
41   PUBLIC   sbc_blk_core       ! routine called in sbcmod module
42   PUBLIC   blk_ice_core       ! routine called in sbc_ice_lim module
43     
44   INTEGER , PARAMETER ::   jpfld   = 8           ! maximum number of files to read
45   INTEGER , PARAMETER ::   jp_wndi = 1           ! index of 10m wind velocity (i-component) (m/s)    at T-point
46   INTEGER , PARAMETER ::   jp_wndj = 2           ! index of 10m wind velocity (j-component) (m/s)    at T-point
47   INTEGER , PARAMETER ::   jp_humi = 3           ! index of specific humidity               ( % )
48   INTEGER , PARAMETER ::   jp_qsr  = 4           ! index of solar heat                      (W/m2)
49   INTEGER , PARAMETER ::   jp_qlw  = 5           ! index of Long wave                       (W/m2)
50   INTEGER , PARAMETER ::   jp_tair = 6           ! index of 10m air temperature             (Kelvin)
51   INTEGER , PARAMETER ::   jp_prec = 7           ! index of total precipitation (rain+snow) (Kg/m2/s)
52   INTEGER , PARAMETER ::   jp_snow = 8           ! index of snow (solid prcipitation)       (kg/m2/s)
53   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf   ! structure of input fields (file informations, fields read)
54         
55   !! * CORE bulk parameters
56   REAL(wp), PARAMETER ::   rhoa =    1.22        ! air density
57   REAL(wp), PARAMETER ::   cpa  = 1000.5         ! specific heat of air
58   REAL(wp), PARAMETER ::   Lv   =    2.5e6       ! latent heat of vaporization
59   REAL(wp), PARAMETER ::   Ls   =    2.839e6     ! latent heat of sublimation
60   REAL(wp), PARAMETER ::   Stef =    5.67e-8     ! Stefan Boltzmann constant
61   REAL(wp), PARAMETER ::   Cice =    1.63e-3     ! transfer coefficient over ice
62
63   LOGICAL  ::   ln_2m = .FALSE.                  !: logical flag for height of air temp. and hum
64   REAL(wp) ::   alpha_precip=1.                  !: multiplication factor for precipitation
65
66   !! * Substitutions
67#  include "domzgr_substitute.h90"
68#  include "vectopt_loop_substitute.h90"
69   !!----------------------------------------------------------------------
70   !!   OPA 9.0 , LOCEAN-IPSL (2006)
71   !! $ Id: $
72   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
73   !!----------------------------------------------------------------------
74
75CONTAINS
76
77   SUBROUTINE sbc_blk_core( kt )
78      !!---------------------------------------------------------------------
79      !!                    ***  ROUTINE sbc_blk_core  ***
80      !!                   
81      !! ** Purpose :   provide at each time step the surface ocean fluxes
82      !!      (momentum, heat, freshwater and runoff)
83      !!
84      !! ** Method  :   READ each fluxes in NetCDF files
85      !!      The i-component of the stress                utau   (N/m2)
86      !!      The j-component of the stress                vtau   (N/m2)
87      !!      the net downward heat flux                   qtot   (watt/m2)
88      !!      the net downward radiative flux              qsr    (watt/m2)
89      !!      the net upward water (evapo - precip)        emp    (kg/m2/s)
90      !!                Assumptions made:
91      !!       - each file content an entire year (read record, not the time axis)
92      !!       - first and last record are part of the previous and next year
93      !!         (useful for time interpolation)
94      !!       - the number of records is 2 + 365*24 / freqh(jf)
95      !!         or 366 in leap year case
96      !!
97      !!      C A U T I O N : never mask the surface stress fields
98      !!                      the stress is assumed to be in the mesh referential
99      !!                      i.e. the (i,j) referential
100      !!
101      !! ** Action  :   defined at each time-step at the air-sea interface
102      !!              - utau  &  vtau   : stress components in geographical ref.
103      !!              - qns   &  qsr    : non solar and solar heat fluxes
104      !!              - emp             : evap - precip (volume flux)
105      !!              - emps            : evap - precip (concentration/dillution)
106      !!----------------------------------------------------------------------
107      INTEGER, INTENT( in  ) ::   kt   ! ocean time step
108      !!
109      INTEGER  ::   jf       ! dummy indices
110      INTEGER  ::   ierror   ! return error code
111      !!
112      CHARACTER(len=100) ::  cn_dir   !   Root directory for location of core files
113      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i     ! array of namelist informations on the fields to read
114      TYPE(FLD_N) ::   sn_wndi, sn_wndj, sn_humi, sn_qsr       ! informations about the fields to be read
115      TYPE(FLD_N) ::   sn_qlw , sn_tair, sn_prec, sn_snow      !   "                                 "
116      NAMELIST/namsbc_core/ cn_dir, ln_2m, alpha_precip, sn_wndi, sn_wndj, sn_humi, sn_qsr,   &
117         &                                               sn_qlw , sn_tair, sn_prec, sn_snow
118      !!---------------------------------------------------------------------
119
120      !                                         ! ====================== !
121      IF( kt == nit000 ) THEN                   !  First call kt=nit000  !
122         !                                      ! ====================== !
123         ! set file information (default values)
124         cn_dir = './'       ! directory in which the model is executed
125
126         ! (NB: frequency positive => hours, negative => months)
127         !            !    file     ! frequency !  variable  ! time intep !  clim  ! starting !
128         !            !    name     !  (hours)  !   name     !   (T/F)    !  (0/1) !  record  !
129         sn_wndi = FLD_N( 'uwnd10m' ,    24.    ,  'u_10'    ,  .FALSE.   ,    0   ,     0    ) 
130         sn_wndj = FLD_N( 'vwnd10m' ,    24.    ,  'v_10'    ,  .FALSE.   ,    0   ,     0    ) 
131         sn_qsr  = FLD_N( 'qsw'     ,    24.    ,  'qsw'     ,  .FALSE.   ,    0   ,     0    ) 
132         sn_qlw  = FLD_N( 'qlw'     ,    24.    ,  'qlw'     ,  .FALSE.   ,    0   ,     0    ) 
133         sn_tair = FLD_N( 'tair10m' ,    24.    ,  't_10'    ,  .FALSE.   ,    0   ,     0    ) 
134         sn_humi = FLD_N( 'humi10m' ,    24.    ,  'q_10'    ,  .FALSE.   ,    0   ,     0    ) 
135         sn_prec = FLD_N( 'precip'  ,   -12.    ,  'precip'  ,  .TRUE.    ,    0   ,     0    ) 
136         sn_snow = FLD_N( 'snow'    ,   -12.    ,  'snow'    ,  .TRUE.    ,    0   ,     0    ) 
137
138         REWIND( numnam )                    ! ... read in namlist namsbc_core
139         READ  ( numnam, namsbc_core )
140
141         ! store namelist information in an array
142         slf_i(jp_wndi) = sn_wndi   ;   slf_i(jp_wndj) = sn_wndj
143         slf_i(jp_qsr ) = sn_qsr    ;   slf_i(jp_qlw ) = sn_qlw
144         slf_i(jp_tair) = sn_tair   ;   slf_i(jp_humi) = sn_humi
145         slf_i(jp_prec) = sn_prec   ;   slf_i(jp_snow) = sn_snow
146         
147         ! set sf structure
148         ALLOCATE( sf(jpfld), STAT=ierror )
149         IF( ierror > 0 ) THEN
150            CALL ctl_stop( 'sbc_blk_core: unable to allocate sf structure' )   ;   RETURN
151         ENDIF
152
153         DO jf = 1, jpfld
154            WRITE(sf(jf)%clrootname,'(a,a)' )   TRIM( cn_dir ), TRIM( slf_i(jf)%clname )
155            sf(jf)%freqh   = slf_i(jf)%freqh
156            sf(jf)%clvar   = slf_i(jf)%clvar
157            sf(jf)%ln_tint = slf_i(jf)%ln_tint
158            sf(jf)%nclim   = slf_i(jf)%nclim
159            sf(jf)%nstrec  = slf_i(jf)%nstrec
160         END DO
161
162         IF(lwp) THEN      ! control print
163            WRITE(numout,*)           
164            WRITE(numout,*) 'sbc_blk_core : flux formulattion for ocean surface boundary condition'
165            WRITE(numout,*) '~~~~~~~~~~~~ '
166            WRITE(numout,*) '          namsbc_core Namelist'
167            WRITE(numout,*) '          ln_2m        = ', ln_2m
168            WRITE(numout,*) '          alpha_precip = ', alpha_precip
169            WRITE(numout,*) '          list of files and frequency (>0: in hours ; <0 in months)'
170            DO jf = 1, jpfld
171                WRITE(numout,*) '               file root name: ' , TRIM( sf(jf)%clrootname ),   &
172                   &                          ' variable name: '  , TRIM( sf(jf)%clvar      )
173                WRITE(numout,*) '               frequency: '      ,       sf(jf)%freqh       ,   &
174                   &                          ' time interp: '    ,       sf(jf)%ln_tint     ,   &
175                   &                          ' climatology: '    ,       sf(jf)%nclim       ,   &
176                   &                          ' starting record: ',       sf(jf)%nstrec
177            END DO
178            IF( ln_2m )   THEN
179               WRITE(numout,*) '          Calling TURB_CORE_2Z for bulk transfert coefficients'
180            ELSE
181               WRITE(numout,*) '          Calling TURB_CORE_1Z for bulk transfert coefficients'
182            ENDIF
183            WRITE(numout,*)           
184            !
185         ENDIF
186         !
187      ENDIF
188
189
190      CALL fld_read( kt, nn_fsbc, sf )                ! input fields provided at the current time-step
191
192#if defined key_lim3
193      tatm_ice(:,:) = sf(jp_tair)%fnow(:,:)
194#endif
195
196      IF( MOD( kt - 1, nn_fsbc ) == 0 ) THEN
197          CALL blk_oce_core( sst_m, ssu_m, ssv_m )    ! compute the surface ocean fluxes using CLIO bulk formulea
198      ENDIF
199      !                                               ! using CORE bulk formulea
200   END SUBROUTINE sbc_blk_core
201   
202   
203   SUBROUTINE blk_oce_core( pst, pu, pv )
204      !!---------------------------------------------------------------------
205      !!                     ***  ROUTINE blk_core  ***
206      !!
207      !! ** Purpose :   provide the momentum, heat and freshwater fluxes at
208      !!      the ocean surface at each time step
209      !!
210      !! ** Method  :   CORE bulk formulea for the ocean using atmospheric
211      !!      fields read in sbc_read
212      !!
213      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
214      !!              - vtau    : j-component of the stress at V-point  (N/m2)
215      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
216      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
217      !!              - evap    : Evaporation over the ocean            (kg/m2/s)
218      !!              - tprecip : Total precipitation                   (Kg/m2/s)
219      !!              - sprecip : Solid precipitation                   (Kg/m2/s)
220      !!---------------------------------------------------------------------
221      REAL(wp), INTENT(in), DIMENSION(jpi,jpj) ::   pst   ! surface temperature                      [Celcius]
222      REAL(wp), INTENT(in), DIMENSION(jpi,jpj) ::   pu    ! surface current at U-point (i-component) [m/s]
223      REAL(wp), INTENT(in), DIMENSION(jpi,jpj) ::   pv    ! surface current at V-point (j-component) [m/s]
224
225      INTEGER  ::   ji, jj     ! dummy loop indices
226      REAL(wp) ::   zcoef_qsatw
227      REAL(wp), DIMENSION(jpi,jpj) ::   zwnd_i, zwnd_j    ! wind speed components at T-point
228      REAL(wp), DIMENSION(jpi,jpj) ::   zqsatw            ! specific humidity at pst
229      REAL(wp), DIMENSION(jpi,jpj) ::   zqlw, zqsb        ! long wave and sensible heat fluxes
230      REAL(wp), DIMENSION(jpi,jpj) ::   zqla, zevap       ! latent heat fluxes and evaporation
231      REAL(wp), DIMENSION(jpi,jpj) ::   Cd                ! transfer coefficient for momentum      (tau)
232      REAL(wp), DIMENSION(jpi,jpj) ::   Ch                ! transfer coefficient for sensible heat (Q_sens)
233      REAL(wp), DIMENSION(jpi,jpj) ::   Ce                ! tansfert coefficient for evaporation   (Q_lat)
234      REAL(wp), DIMENSION(jpi,jpj) ::   zst               ! surface temperature in Kelvin
235      REAL(wp), DIMENSION(jpi,jpj) ::   zt_zu             ! air temperature at wind speed height
236      REAL(wp), DIMENSION(jpi,jpj) ::   zq_zu             ! air spec. hum.  at wind speed height
237      !!---------------------------------------------------------------------
238
239      ! local scalars ( place there for vector optimisation purposes)
240      zcoef_qsatw = 0.98 * 640380. / rhoa
241     
242      zst(:,:) = pst(:,:) + rt0      ! converte Celcius to Kelvin (and set minimum value far above 0 K)
243
244      ! ----------------------------------------------------------------------------- !
245      !      0   Wind components and module at T-point relative to the moving ocean   !
246      ! ----------------------------------------------------------------------------- !
247
248      ! ... components ( U10m - U_oce ) at T-point (unmasked)
249      zwnd_i(:,:) = 0.e0 
250      zwnd_j(:,:) = 0.e0
251#if defined key_vectopt_loop
252!CDIR COLLAPSE
253#endif
254      DO jj = 2, jpjm1
255         DO ji = fs_2, fs_jpim1   ! vect. opt.
256            zwnd_i(ji,jj) = (  sf(jp_wndi)%fnow(ji,jj) - 0.5 * ( pu(ji-1,jj  ) + pu(ji,jj) )  )
257            zwnd_j(ji,jj) = (  sf(jp_wndj)%fnow(ji,jj) - 0.5 * ( pv(ji  ,jj-1) + pv(ji,jj) )  )
258         END DO
259      END DO
260      CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. )
261      CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. )
262      ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked)
263!CDIR NOVERRCHK
264!CDIR COLLAPSE
265      wndm(:,:) = SQRT(  zwnd_i(:,:) * zwnd_i(:,:)   &
266         &             + zwnd_j(:,:) * zwnd_j(:,:)  ) * tmask(:,:,1)
267
268      ! ----------------------------------------------------------------------------- !
269      !      I   Radiative FLUXES                                                     !
270      ! ----------------------------------------------------------------------------- !
271   
272      ! ocean albedo assumed to be 0.066
273!CDIR COLLAPSE
274      qsr (:,:) = ( 1. - 0.066 ) * sf(jp_qsr)%fnow(:,:) * tmask(:,:,1)                                 ! Short Wave
275!CDIR COLLAPSE
276      zqlw(:,:) = (  sf(jp_qlw)%fnow(:,:) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:)  ) * tmask(:,:,1)   ! Long  Wave
277                     
278      ! ----------------------------------------------------------------------------- !
279      !     II    Turbulent FLUXES                                                    !
280      ! ----------------------------------------------------------------------------- !
281
282      ! ... specific humidity at SST and IST
283!CDIR NOVERRCHK
284!CDIR COLLAPSE
285      zqsatw(:,:) = zcoef_qsatw * EXP( -5107.4 / zst(:,:) ) 
286
287      ! ... NCAR Bulk formulae, computation of Cd, Ch, Ce at T-point :
288      IF( ln_2m ) THEN
289         !! If air temp. and spec. hum. are given at different height (2m) than wind (10m) :
290         CALL TURB_CORE_2Z(2.,10., zst   , sf(jp_tair)%fnow,         &
291            &                      zqsatw, sf(jp_humi)%fnow, wndm,   &
292            &                      Cd    , Ch              , Ce  ,   &
293            &                      zt_zu , zq_zu                   )
294      ELSE
295         !! If air temp. and spec. hum. are given at same height than wind (10m) :
296!gm bug?  at the compiling phase, add a copy in temporary arrays...  ==> check perf
297!         CALL TURB_CORE_1Z( 10., zst   (:,:), sf(jp_tair)%fnow(:,:),              &
298!            &                    zqsatw(:,:), sf(jp_humi)%fnow(:,:), wndm(:,:),   &
299!            &                    Cd    (:,:),             Ch  (:,:), Ce  (:,:)  )
300!gm bug
301         CALL TURB_CORE_1Z( 10., zst   , sf(jp_tair)%fnow,         &
302            &                    zqsatw, sf(jp_humi)%fnow, wndm,   &
303            &                    Cd    , Ch              , Ce    )
304      ENDIF
305
306      ! ... utau, vtau at U- and V_points, resp.
307      !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines
308      zwnd_i(:,:) = rhoa * wndm(:,:) * Cd(:,:) * zwnd_i(:,:)
309      zwnd_j(:,:) = rhoa * wndm(:,:) * Cd(:,:) * zwnd_j(:,:)
310      DO jj = 1, jpjm1
311         DO ji = 1, fs_jpim1
312            utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj  ) )
313            vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji  ,jj+1) )
314         END DO
315      END DO
316      CALL lbc_lnk( utau(:,:), 'U', -1. )
317      CALL lbc_lnk( vtau(:,:), 'V', -1. )
318
319      !  Turbulent fluxes over ocean
320      ! -----------------------------
321      IF( ln_2m ) THEN
322         ! Values of temp. and hum. adjusted to 10m must be used instead of 2m values
323         zevap(:,:) = MAX( 0.e0, rhoa    *Ce(:,:)*( zqsatw(:,:) - zq_zu(:,:) ) * wndm(:,:) )   ! Evaporation
324         zqsb (:,:) =            rhoa*cpa*Ch(:,:)*( zst   (:,:) - zt_zu(:,:) ) * wndm(:,:)     ! Sensible Heat
325      ELSE
326!CDIR COLLAPSE
327         zevap(:,:) = MAX( 0.e0, rhoa    *Ce(:,:)*( zqsatw(:,:) - sf(jp_humi)%fnow(:,:) ) * wndm(:,:) )   ! Evaporation
328!CDIR COLLAPSE
329         zqsb (:,:) =            rhoa*cpa*Ch(:,:)*( zst   (:,:) - sf(jp_tair)%fnow(:,:) ) * wndm(:,:)     ! Sensible Heat
330      ENDIF
331!CDIR COLLAPSE
332      zqla (:,:) = Lv * zevap(:,:)                                                              ! Latent Heat
333
334      IF(ln_ctl) THEN
335         CALL prt_ctl( tab2d_1=zqla  , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=Ce , clinfo2=' Ce  : ' )
336         CALL prt_ctl( tab2d_1=zqsb  , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=Ch , clinfo2=' Ch  : ' )
337         CALL prt_ctl( tab2d_1=zqlw  , clinfo1=' blk_oce_core: zqlw   : ', tab2d_2=qsr, clinfo2=' qsr : ' )
338         CALL prt_ctl( tab2d_1=zqsatw, clinfo1=' blk_oce_core: zqsatw : ', tab2d_2=zst, clinfo2=' zst : ' )
339         CALL prt_ctl( tab2d_1=utau  , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
340            &          tab2d_2=vtau  , clinfo2=              ' vtau : '  , mask2=vmask )
341         CALL prt_ctl( tab2d_1=wndm  , clinfo1=' blk_oce_core: wndm   : ')
342         CALL prt_ctl( tab2d_1=zst   , clinfo1=' blk_oce_core: zst    : ')
343      ENDIF
344       
345      ! ----------------------------------------------------------------------------- !
346      !     III    Total FLUXES                                                       !
347      ! ----------------------------------------------------------------------------- !
348     
349!CDIR COLLAPSE
350      qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)      ! Downward Non Solar flux
351!CDIR COLLAPSE
352      emp (:,:) = zevap(:,:) - sf(jp_prec)%fnow(:,:) * alpha_precip * tmask(:,:,1)
353!CDIR COLLAPSE
354      emps(:,:) = zevap(:,:) - sf(jp_prec)%fnow(:,:) * alpha_precip * tmask(:,:,1)
355      !
356   END SUBROUTINE blk_oce_core
357   
358   
359   SUBROUTINE blk_ice_core(  pst   , pui   , pvi   , palb ,   &
360      &                      p_taui, p_tauj, p_qns , p_qsr,   &
361      &                      p_qla , p_dqns, p_dqla,          &
362      &                      p_tpr , p_spr ,                  &
363      &                      p_fr1 , p_fr2 , cd_grid  ) 
364      !!---------------------------------------------------------------------
365      !!                     ***  ROUTINE blk_ice_core  ***
366      !!
367      !! ** Purpose :   provide the surface boundary condition over sea-ice
368      !!
369      !! ** Method  :   compute momentum, heat and freshwater exchanged
370      !!                between atmosphere and sea-ice using CORE bulk
371      !!                formulea, ice variables and read atmmospheric fields.
372      !!                NB: ice drag coefficient is assumed to be a constant
373      !!
374      !! caution : the net upward water flux has with mm/day unit
375      !!---------------------------------------------------------------------
376      REAL(wp), INTENT(in   ), DIMENSION(:,:,:)   ::   pst      ! ice surface temperature (>0, =rt0 over land) [Kelvin]
377      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   pui      ! ice surface velocity (i- and i- components      [m/s]
378      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   pvi      !    at I-point (B-grid) or U & V-point (C-grid)
379      REAL(wp), INTENT(in   ), DIMENSION(:,:,:)   ::   palb     ! ice albedo (clear sky) (alb_ice_cs)               [%]
380      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   p_taui   ! i- & j-components of surface ice stress        [N/m2]
381      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   p_tauj   !   at I-point (B-grid) or U & V-point (C-grid)
382      REAL(wp), INTENT(  out), DIMENSION(:,:,:)   ::   p_qns    ! non solar heat flux over ice (T-point)         [W/m2]
383      REAL(wp), INTENT(  out), DIMENSION(:,:,:)   ::   p_qsr    !     solar heat flux over ice (T-point)         [W/m2]
384      REAL(wp), INTENT(  out), DIMENSION(:,:,:)   ::   p_qla    ! latent    heat flux over ice (T-point)         [W/m2]
385      REAL(wp), INTENT(  out), DIMENSION(:,:,:)   ::   p_dqns   ! non solar heat sensistivity  (T-point)         [W/m2]
386      REAL(wp), INTENT(  out), DIMENSION(:,:,:)   ::   p_dqla   ! latent    heat sensistivity  (T-point)         [W/m2]
387      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   p_tpr    ! total precipitation          (T-point)      [Kg/m2/s]
388      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   p_spr    ! solid precipitation          (T-point)      [Kg/m2/s]
389      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   p_fr1    ! 1sr fraction of qsr penetration in ice (T-point)  [%]
390      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   p_fr2    ! 2nd fraction of qsr penetration in ice (T-point)  [%]
391      CHARACTER(len=1), INTENT(in   )             ::   cd_grid  ! ice grid ( C or B-grid)
392      !!
393      INTEGER  ::   ji, jj, jl    ! dummy loop indices
394      INTEGER  ::   ijpl          ! number of ice categories (size of 3rd dim of input arrays)
395      REAL(wp) ::   zst2, zst3
396      REAL(wp) ::   zcoef_wnorm, zcoef_wnorm2, zcoef_dqlw, zcoef_dqla, zcoef_dqsb
397      REAL(wp) ::   zcoef_frca                                   ! fractional cloud amount
398      REAL(wp) ::   zwnorm_f, zwndi_f , zwndj_f                  ! relative wind module and components at F-point
399      REAL(wp) ::             zwndi_t , zwndj_t                  ! relative wind components at T-point
400      REAL(wp), DIMENSION(jpi,jpj) ::   z_wnds_t                 ! wind speed ( = | U10m - U_ice | ) at T-point
401      REAL(wp), DIMENSION(jpi,jpj,SIZE(pst,3)) ::   z_qlw        ! long wave heat flux over ice
402      REAL(wp), DIMENSION(jpi,jpj,SIZE(pst,3)) ::   z_qsb        ! sensible  heat flux over ice
403      REAL(wp), DIMENSION(jpi,jpj,SIZE(pst,3)) ::   z_dqlw       ! sensible  heat flux over ice
404      REAL(wp), DIMENSION(jpi,jpj,SIZE(pst,3)) ::   z_dqsb       ! sensible  heat flux over ice
405      !!---------------------------------------------------------------------
406
407      ijpl  = SIZE( pst, 3 )                 ! number of ice categories
408
409      ! local scalars ( place there for vector optimisation purposes)
410      zcoef_wnorm = rhoa * Cice
411      zcoef_wnorm2 = rhoa * Cice * 0.5
412      zcoef_dqlw = 4.0 * 0.95 * Stef
413      zcoef_dqla = -Ls * Cice * 11637800. * (-5897.8)
414      zcoef_dqsb = rhoa * cpa * Cice
415      zcoef_frca = 1.0  - 0.3
416
417!!gm brutal....
418      z_wnds_t(:,:) = 0.e0
419      p_taui  (:,:) = 0.e0
420      p_tauj  (:,:) = 0.e0
421!!gm end
422
423      ! ----------------------------------------------------------------------------- !
424      !    Wind components and module relative to the moving ocean ( U10m - U_ice )   !
425      ! ----------------------------------------------------------------------------- !
426      SELECT CASE( cd_grid )
427      CASE( 'B' )                  ! B-grid ice dynamics :   I-point (i.e. F-point with sea-ice indexation)
428         !                           and scalar wind at T-point ( = | U10m - U_ice | ) (masked)
429#if defined key_vectopt_loop
430!CDIR COLLAPSE
431#endif
432!CDIR NOVERRCHK
433         DO jj = 2, jpjm1
434            DO ji = fs_2, fs_jpim1
435               ! ... scalar wind at I-point (fld being at T-point)
436               zwndi_f = 0.25 * (  sf(jp_wndi)%fnow(ji-1,jj  ) + sf(jp_wndi)%fnow(ji  ,jj  )   &
437                  &              + sf(jp_wndi)%fnow(ji-1,jj-1) + sf(jp_wndi)%fnow(ji  ,jj-1)  ) - pui(ji,jj)
438               zwndj_f = 0.25 * (  sf(jp_wndj)%fnow(ji-1,jj  ) + sf(jp_wndj)%fnow(ji  ,jj  )   &
439                  &              + sf(jp_wndj)%fnow(ji-1,jj-1) + sf(jp_wndj)%fnow(ji  ,jj-1)  ) - pvi(ji,jj)
440               zwnorm_f = zcoef_wnorm * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f )
441               ! ... ice stress at I-point
442               p_taui(ji,jj) = zwnorm_f * zwndi_f
443               p_tauj(ji,jj) = zwnorm_f * zwndj_f
444               ! ... scalar wind at T-point (fld being at T-point)
445               zwndi_t = sf(jp_wndi)%fnow(ji,jj) - 0.25 * (  pui(ji,jj+1) + pui(ji+1,jj+1)   &
446                  &                                        + pui(ji,jj  ) + pui(ji+1,jj  )  )
447               zwndj_t = sf(jp_wndj)%fnow(ji,jj) - 0.25 * (  pvi(ji,jj+1) + pvi(ji+1,jj+1)   &
448                  &                                        + pvi(ji,jj  ) + pvi(ji+1,jj  )  )
449               z_wnds_t(ji,jj)  = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
450            END DO
451         END DO
452         CALL lbc_lnk( p_taui  , 'I', -1. )
453         CALL lbc_lnk( p_tauj  , 'I', -1. )
454         CALL lbc_lnk( z_wnds_t, 'T',  1. )
455         !
456      CASE( 'C' )                  ! C-grid ice dynamics :   U & V-points (same as ocean)
457#if defined key_vectopt_loop
458!CDIR COLLAPSE
459#endif
460         DO jj = 2, jpj
461            DO ji = fs_2, jpi   ! vect. opt.
462               zwndi_t = (  sf(jp_wndi)%fnow(ji,jj) - 0.5 * ( pui(ji-1,jj  ) + pui(ji,jj) )  )
463               zwndj_t = (  sf(jp_wndj)%fnow(ji,jj) - 0.5 * ( pvi(ji  ,jj-1) + pvi(ji,jj) )  )
464               z_wnds_t(ji,jj)  = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
465            END DO
466         END DO
467#if defined key_vectopt_loop
468!CDIR COLLAPSE
469#endif
470         DO jj = 2, jpjm1
471            DO ji = fs_2, fs_jpim1   ! vect. opt.
472               p_taui(ji,jj) = zcoef_wnorm2 * ( z_wnds_t(ji+1,jj) + z_wnds_t(ji,jj) )                          &
473                  &          * ( 0.5 * (sf(jp_wndi)%fnow(ji+1,jj) + sf(jp_wndi)%fnow(ji,jj) ) - pui(ji,jj) )
474               p_tauj(ji,jj) = zcoef_wnorm2 * ( z_wnds_t(ji,jj+1) + z_wnds_t(ji,jj) )                          &
475                  &          * ( 0.5 * (sf(jp_wndj)%fnow(ji,jj+1) + sf(jp_wndj)%fnow(ji,jj) ) - pvi(ji,jj) )
476            END DO
477         END DO
478         CALL lbc_lnk( p_taui  , 'U', -1. )
479         CALL lbc_lnk( p_tauj  , 'V', -1. )
480         CALL lbc_lnk( z_wnds_t, 'T',  1. )
481         !
482      END SELECT
483
484      !                                     ! ========================== !
485      DO jl = 1, ijpl                       !  Loop over ice categories  !
486         !                                  ! ========================== !
487!CDIR NOVERRCHK
488!CDIR COLLAPSE
489         DO jj = 1 , jpj
490!CDIR NOVERRCHK
491            DO ji = 1, jpi
492               ! ----------------------------!
493               !      I   Radiative FLUXES   !
494               ! ----------------------------!
495               zst2 = pst(ji,jj,jl) * pst(ji,jj,jl)
496               zst3 = pst(ji,jj,jl) * zst2
497               ! Short Wave (sw)
498               p_qsr(ji,jj,jl) = ( 1. - palb(ji,jj,jl) ) * sf(jp_qsr)%fnow(ji,jj) * tmask(ji,jj,1)
499               ! Long  Wave (lw)
500               z_qlw(ji,jj,jl) = 0.95 * (  sf(jp_qlw)%fnow(ji,jj)       &                         
501                  &                   - Stef * pst(ji,jj,jl) * zst3  ) * tmask(ji,jj,1)
502               ! lw sensitivity
503               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3                                               
504
505               ! ----------------------------!
506               !     II    Turbulent FLUXES  !
507               ! ----------------------------!
508
509               ! ... turbulent heat fluxes
510               ! Sensible Heat
511               z_qsb(ji,jj,jl) = rhoa * cpa * Cice * z_wnds_t(ji,jj) * ( pst(ji,jj,jl) - sf(jp_tair)%fnow(ji,jj) )
512               ! Latent Heat
513               p_qla(ji,jj,jl) = MAX( 0.e0, rhoa * Ls  * Cice * z_wnds_t(ji,jj)   &                           
514                  &                    * (  11637800. * EXP( -5897.8 / pst(ji,jj,jl) ) / rhoa - sf(jp_humi)%fnow(ji,jj)  ) )
515               ! Latent heat sensitivity for ice (Dqla/Dt)
516               p_dqla(ji,jj,jl) = zcoef_dqla * z_wnds_t(ji,jj) / ( zst2 ) * EXP( -5897.8 / pst(ji,jj,jl) )
517               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
518               z_dqsb(ji,jj,jl) = zcoef_dqsb * z_wnds_t(ji,jj)
519
520               ! ----------------------------!
521               !     III    Total FLUXES     !
522               ! ----------------------------!
523               ! Downward Non Solar flux
524               p_qns (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - p_qla (ji,jj,jl)     
525               ! Total non solar heat flux sensitivity for ice
526               p_dqns(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + p_dqla(ji,jj,jl) )   
527            END DO
528            !
529         END DO
530         !
531      END DO
532      !
533      !--------------------------------------------------------------------
534      ! FRACTIONs of net shortwave radiation which is not absorbed in the
535      ! thin surface layer and penetrates inside the ice cover
536      ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 )
537   
538!CDIR COLLAPSE
539      p_fr1(:,:) = ( 0.18 * ( 1.0 - zcoef_frca ) + 0.35 * zcoef_frca )
540!CDIR COLLAPSE
541      p_fr2(:,:) = ( 0.82 * ( 1.0 - zcoef_frca ) + 0.65 * zcoef_frca )
542       
543!CDIR COLLAPSE
544      p_tpr(:,:) = sf(jp_prec)%fnow(:,:) * alpha_precip      ! total precipitation [kg/m2/s]
545!CDIR COLLAPSE
546      p_spr(:,:) = sf(jp_snow)%fnow(:,:) * alpha_precip      ! solid precipitation [kg/m2/s]
547      !
548      IF(ln_ctl) THEN
549         CALL prt_ctl(tab3d_1=p_qla   , clinfo1=' blk_ice_core: p_qla  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb  : ', kdim=ijpl)
550         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice_core: z_qlw  : ', tab3d_2=p_dqla  , clinfo2=' p_dqla : ', kdim=ijpl)
551         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice_core: z_dqsb : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw : ', kdim=ijpl)
552         CALL prt_ctl(tab3d_1=p_dqns  , clinfo1=' blk_ice_core: p_dqns : ', tab3d_2=p_qsr   , clinfo2=' p_qsr  : ', kdim=ijpl)
553         CALL prt_ctl(tab3d_1=pst     , clinfo1=' blk_ice_core: pst    : ', tab3d_2=p_qns   , clinfo2=' p_qns  : ', kdim=ijpl)
554         CALL prt_ctl(tab2d_1=p_tpr   , clinfo1=' blk_ice_core: p_tpr  : ', tab2d_2=p_spr   , clinfo2=' p_spr  : ')
555         CALL prt_ctl(tab2d_1=p_taui  , clinfo1=' blk_ice_core: p_taui : ', tab2d_2=p_tauj  , clinfo2=' p_tauj : ')
556         CALL prt_ctl(tab2d_1=z_wnds_t, clinfo1=' blk_ice_core: z_wnds_t : ')
557      ENDIF
558
559   END SUBROUTINE blk_ice_core
560 
561
562   SUBROUTINE TURB_CORE_1Z(zu, sst, T_a, q_sat, q_a,   &
563      &                        dU, Cd, Ch, Ce   )
564      !!----------------------------------------------------------------------
565      !!                      ***  ROUTINE  turb_core  ***
566      !!
567      !! ** Purpose :   Computes turbulent transfert coefficients of surface
568      !!                fluxes according to Large & Yeager (2004)
569      !!
570      !! ** Method  :   I N E R T I A L   D I S S I P A T I O N   M E T H O D
571      !!      Momentum, Latent and sensible heat exchange coefficients
572      !!      Caution: this procedure should only be used in cases when air
573      !!      temperature (T_air), air specific humidity (q_air) and wind (dU)
574      !!      are provided at the same height 'zzu'!
575      !!
576      !! References :
577      !!      Large & Yeager, 2004 : ???
578      !! History :
579      !!        !  XX-XX  (???     )  Original code
580      !!   9.0  !  05-08  (L. Brodeau) Rewriting and optimization
581      !!----------------------------------------------------------------------
582      !! * Arguments
583
584      REAL(wp), INTENT(in) :: zu                 ! altitude of wind measurement       [m]
585      REAL(wp), INTENT(in), DIMENSION(jpi,jpj) ::  &
586         sst,       &       ! sea surface temperature         [Kelvin]
587         T_a,       &       ! potential air temperature       [Kelvin]
588         q_sat,     &       ! sea surface specific humidity   [kg/kg]
589         q_a,       &       ! specific air humidity           [kg/kg]
590         dU                 ! wind module |U(zu)-U(0)|        [m/s]
591      REAL(wp), intent(out), DIMENSION(jpi,jpj)  :: &
592         Cd,    &                ! transfert coefficient for momentum       (tau)
593         Ch,    &                ! transfert coefficient for temperature (Q_sens)
594         Ce                      ! transfert coefficient for evaporation  (Q_lat)
595
596      !! * Local declarations
597      REAL(wp), DIMENSION(jpi,jpj)  ::   &
598         dU10,        &       ! dU                                   [m/s]
599         dT,          &       ! air/sea temperature differeence      [K]
600         dq,          &       ! air/sea humidity difference          [K]
601         Cd_n10,      &       ! 10m neutral drag coefficient
602         Ce_n10,      &       ! 10m neutral latent coefficient
603         Ch_n10,      &       ! 10m neutral sensible coefficient
604         sqrt_Cd_n10, &       ! root square of Cd_n10
605         sqrt_Cd,     &       ! root square of Cd
606         T_vpot,      &       ! virtual potential temperature        [K]
607         T_star,      &       ! turbulent scale of tem. fluct.
608         q_star,      &       ! turbulent humidity of temp. fluct.
609         U_star,      &       ! turb. scale of velocity fluct.
610         L,           &       ! Monin-Obukov length                  [m]
611         zeta,        &       ! stability parameter at height zu
612         U_n10,       &       ! neutral wind velocity at 10m         [m]   
613         xlogt, xct, zpsi_h, zpsi_m
614      !!
615      INTEGER :: j_itt
616      INTEGER, PARAMETER :: nb_itt = 3
617      INTEGER, DIMENSION(jpi,jpj)  ::   &
618         stab         ! 1st guess stability test integer
619
620      REAL(wp), PARAMETER ::                        &
621         grav   = 9.8,          &  ! gravity                       
622         kappa  = 0.4              ! von Karman s constant
623
624      !! * Start
625      !! Air/sea differences
626      dU10 = max(0.5, dU)     ! we don't want to fall under 0.5 m/s
627      dT = T_a - sst          ! assuming that T_a is allready the potential temp. at zzu
628      dq = q_a - q_sat
629      !!   
630      !! Virtual potential temperature
631      T_vpot = T_a*(1. + 0.608*q_a)
632      !!
633      !! Neutral Drag Coefficient
634      stab    = 0.5 + sign(0.5,dT)    ! stable : stab = 1 ; unstable : stab = 0
635      Cd_n10  = 1E-3 * ( 2.7/dU10 + 0.142 + dU10/13.09 )    !   L & Y eq. (6a)
636      sqrt_Cd_n10 = sqrt(Cd_n10)
637      Ce_n10  = 1E-3 * ( 34.6 * sqrt_Cd_n10 )               !   L & Y eq. (6b)
638      Ch_n10  = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1-stab)) !   L & Y eq. (6c), (6d)
639      !!
640      !! Initializing transfert coefficients with their first guess neutral equivalents :
641      Cd = Cd_n10 ;  Ce = Ce_n10 ;  Ch = Ch_n10 ;  sqrt_Cd = sqrt(Cd)
642
643      !! * Now starting iteration loop
644      DO j_itt=1, nb_itt
645         !! Turbulent scales :
646         U_star  = sqrt_Cd*dU10                !   L & Y eq. (7a)
647         T_star  = Ch/sqrt_Cd*dT               !   L & Y eq. (7b)
648         q_star  = Ce/sqrt_Cd*dq               !   L & Y eq. (7c)
649
650         !! Estimate the Monin-Obukov length :
651         L  = (U_star**2)/( kappa*grav*(T_star/T_vpot + q_star/(q_a + 1./0.608)) )
652
653         !! Stability parameters :
654         zeta = zu/L ;  zeta   = sign( min(abs(zeta),10.0), zeta )
655         zpsi_h = psi_h(zeta)
656         zpsi_m = psi_m(zeta)
657
658         !! Shifting the wind speed to 10m and neutral stability :
659         U_n10 = dU10*1./(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)) !  L & Y eq. (9a)
660
661         !! Updating the neutral 10m transfer coefficients :
662         Cd_n10  = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09)              !  L & Y eq. (6a)
663         sqrt_Cd_n10 = sqrt(Cd_n10)
664         Ce_n10  = 1E-3 * (34.6 * sqrt_Cd_n10)                           !  L & Y eq. (6b)
665         stab    = 0.5 + sign(0.5,zeta)
666         Ch_n10  = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab))           !  L & Y eq. (6c), (6d)
667
668         !! Shifting the neutral  10m transfer coefficients to ( zu , zeta ) :
669         !!
670         xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10) - zpsi_m)
671         Cd  = Cd_n10/(xct*xct) ;  sqrt_Cd = sqrt(Cd)
672         !!
673         xlogt = log(zu/10.) - zpsi_h
674         !!
675         xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10
676         Ch  = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct
677         !!
678         xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10
679         Ce  = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct
680         !!
681      END DO
682      !!
683    END SUBROUTINE TURB_CORE_1Z
684
685
686    SUBROUTINE TURB_CORE_2Z(zt, zu, sst, T_zt, q_sat, q_zt, dU, Cd, Ch, Ce, T_zu, q_zu)
687      !!----------------------------------------------------------------------
688      !!                      ***  ROUTINE  turb_core  ***
689      !!
690      !! ** Purpose :   Computes turbulent transfert coefficients of surface
691      !!                fluxes according to Large & Yeager (2004).
692      !!
693      !! ** Method  :   I N E R T I A L   D I S S I P A T I O N   M E T H O D
694      !!      Momentum, Latent and sensible heat exchange coefficients
695      !!      Caution: this procedure should only be used in cases when air
696      !!      temperature (T_air) and air specific humidity (q_air) are at 2m
697      !!      whereas wind (dU) is at 10m.
698      !!
699      !! References :
700      !!      Large & Yeager, 2004 : ???
701      !! History :
702      !!   9.0  !  06-12  (L. Brodeau) Original code for 2Z
703      !!----------------------------------------------------------------------
704      !! * Arguments
705      REAL(wp), INTENT(in)   :: &
706         zt,      &     ! height for T_zt and q_zt                   [m]
707         zu             ! height for dU                              [m]
708      REAL(wp), INTENT(in), DIMENSION(jpi,jpj) ::  &
709         sst,      &     ! sea surface temperature              [Kelvin]
710         T_zt,     &     ! potential air temperature            [Kelvin]
711         q_sat,    &     ! sea surface specific humidity         [kg/kg]
712         q_zt,     &     ! specific air humidity                 [kg/kg]
713         dU              ! relative wind module |U(zu)-U(0)|       [m/s]
714      REAL(wp), INTENT(out), DIMENSION(jpi,jpj)  ::  &
715         Cd,       &     ! transfer coefficient for momentum         (tau)
716         Ch,       &     ! transfer coefficient for sensible heat (Q_sens)
717         Ce,       &     ! transfert coefficient for evaporation   (Q_lat)
718         T_zu,     &     ! air temp. shifted at zu                     [K]
719         q_zu            ! spec. hum.  shifted at zu               [kg/kg]
720
721      !! * Local declarations
722      REAL(wp), DIMENSION(jpi,jpj) ::  &
723         dU10,        &     ! dU                                [m/s]
724         dT,          &     ! air/sea temperature differeence   [K]
725         dq,          &     ! air/sea humidity difference       [K]
726         Cd_n10,      &     ! 10m neutral drag coefficient
727         Ce_n10,      &     ! 10m neutral latent coefficient
728         Ch_n10,      &     ! 10m neutral sensible coefficient
729         sqrt_Cd_n10, &     ! root square of Cd_n10
730         sqrt_Cd,     &     ! root square of Cd
731         T_vpot_u,    &     ! virtual potential temperature        [K]
732         T_star,      &     ! turbulent scale of tem. fluct.
733         q_star,      &     ! turbulent humidity of temp. fluct.
734         U_star,      &     ! turb. scale of velocity fluct.
735         L,           &     ! Monin-Obukov length                  [m]
736         zeta_u,      &     ! stability parameter at height zu
737         zeta_t,      &     ! stability parameter at height zt
738         U_n10,       &     ! neutral wind velocity at 10m        [m]
739         xlogt, xct, zpsi_hu, zpsi_ht, zpsi_m
740
741      INTEGER :: j_itt
742      INTEGER, PARAMETER :: nb_itt = 3   ! number of itterations
743      INTEGER, DIMENSION(jpi,jpj) :: &
744           &     stab                ! 1st stability test integer
745      REAL(wp), PARAMETER ::                        &
746         grav   = 9.8,      &  ! gravity                       
747         kappa  = 0.4          ! von Karman's constant
748
749      !!  * Start
750
751      !! Initial air/sea differences
752      dU10 = max(0.5, dU)      !  we don't want to fall under 0.5 m/s
753      dT = T_zt - sst 
754      dq = q_zt - q_sat
755
756      !! Neutral Drag Coefficient :
757      stab = 0.5 + sign(0.5,dT)                 ! stab = 1  if dT > 0  -> STABLE
758      Cd_n10  = 1E-3*( 2.7/dU10 + 0.142 + dU10/13.09 ) 
759      sqrt_Cd_n10 = sqrt(Cd_n10)
760      Ce_n10  = 1E-3*( 34.6 * sqrt_Cd_n10 )
761      Ch_n10  = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1 - stab))
762
763      !! Initializing transf. coeff. with their first guess neutral equivalents :
764      Cd = Cd_n10 ;  Ce = Ce_n10 ;  Ch = Ch_n10 ;  sqrt_Cd = sqrt(Cd)
765
766      !! Initializing z_u values with z_t values :
767      T_zu = T_zt ;  q_zu = q_zt
768
769      !!  * Now starting iteration loop
770      DO j_itt=1, nb_itt
771         dT = T_zu - sst ;  dq = q_zu - q_sat ! Updating air/sea differences
772         T_vpot_u = T_zu*(1. + 0.608*q_zu)    ! Updating virtual potential temperature at zu
773         U_star = sqrt_Cd*dU10                ! Updating turbulent scales :   (L & Y eq. (7))
774         T_star  = Ch/sqrt_Cd*dT              !
775         q_star  = Ce/sqrt_Cd*dq              !
776         !!
777         L = (U_star*U_star) &                ! Estimate the Monin-Obukov length at height zu
778              & / (kappa*grav/T_vpot_u*(T_star*(1.+0.608*q_zu) + 0.608*T_zu*q_star))
779         !! Stability parameters :
780         zeta_u  = zu/L  ;  zeta_u = sign( min(abs(zeta_u),10.0), zeta_u )
781         zeta_t  = zt/L  ;  zeta_t = sign( min(abs(zeta_t),10.0), zeta_t )
782         zpsi_hu = psi_h(zeta_u)
783         zpsi_ht = psi_h(zeta_t)
784         zpsi_m  = psi_m(zeta_u)
785         !!
786         !! Shifting the wind speed to 10m and neutral stability : (L & Y eq.(9a))
787!        U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u)))
788         U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m))
789         !!
790         !! Shifting temperature and humidity at zu :          (L & Y eq. (9b-9c))
791!        T_zu = T_zt - T_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t))
792         T_zu = T_zt - T_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht)
793!        q_zu = q_zt - q_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t))
794         q_zu = q_zt - q_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht)
795         !!
796         !! q_zu cannot have a negative value : forcing 0
797         stab = 0.5 + sign(0.5,q_zu) ;  q_zu = stab*q_zu
798         !!
799         !! Updating the neutral 10m transfer coefficients :
800         Cd_n10  = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09)    ! L & Y eq. (6a)
801         sqrt_Cd_n10 = sqrt(Cd_n10)
802         Ce_n10  = 1E-3 * (34.6 * sqrt_Cd_n10)                 ! L & Y eq. (6b)
803         stab    = 0.5 + sign(0.5,zeta_u)
804         Ch_n10  = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab)) ! L & Y eq. (6c-6d)
805         !!
806         !!
807         !! Shifting the neutral 10m transfer coefficients to (zu,zeta_u) :
808!        xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u))
809         xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)
810         Cd = Cd_n10/(xct*xct) ; sqrt_Cd = sqrt(Cd)
811         !!
812!        xlogt = log(zu/10.) - psi_h(zeta_u)
813         xlogt = log(zu/10.) - zpsi_hu
814         !!
815         xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10
816         Ch  = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct
817         !!
818         xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10
819         Ce  = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct
820         !!
821         !!
822      END DO
823      !!
824    END SUBROUTINE TURB_CORE_2Z
825
826
827    FUNCTION psi_m(zta)   !! Psis, L & Y eq. (8c), (8d), (8e)
828      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: zta
829
830      REAL(wp), PARAMETER :: pi = 3.141592653589793_wp
831      REAL(wp), DIMENSION(jpi,jpj)             :: psi_m
832      REAL(wp), DIMENSION(jpi,jpj)             :: X2, X, stabit
833      X2 = sqrt(abs(1. - 16.*zta))  ;  X2 = max(X2 , 1.0) ;  X  = sqrt(X2)
834      stabit    = 0.5 + sign(0.5,zta)
835      psi_m = -5.*zta*stabit  &                                                  ! Stable
836           & + (1. - stabit)*(2*log((1. + X)/2) + log((1. + X2)/2) - 2*atan(X) + pi/2)  ! Unstable
837    END FUNCTION psi_m
838
839    FUNCTION psi_h(zta)    !! Psis, L & Y eq. (8c), (8d), (8e)
840      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: zta
841
842      REAL(wp), DIMENSION(jpi,jpj)             :: psi_h
843      REAL(wp), DIMENSION(jpi,jpj)             :: X2, X, stabit
844      X2 = sqrt(abs(1. - 16.*zta))  ;  X2 = max(X2 , 1.) ;  X  = sqrt(X2)
845      stabit    = 0.5 + sign(0.5,zta)
846      psi_h = -5.*zta*stabit  &                                       ! Stable
847           & + (1. - stabit)*(2.*log( (1. + X2)/2. ))                 ! Unstable
848    END FUNCTION psi_h
849 
850   !!======================================================================
851END MODULE sbcblk_core
Note: See TracBrowser for help on using the repository browser.