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 @ 888

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

merge dev_001_SBC branche with the trunk to include the New Surface Module package, see ticket: #113

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