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

Last change on this file since 1026 was 1025, checked in by cetlod, 16 years ago

adding wind speed module variable, see ticket 172

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