source: trunk/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 86

Last change on this file since 86 was 86, checked in by cholod, 10 years ago

Print namsbc_core parameters (rn_usecrt...) in ocean.output

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