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 branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 3772

Last change on this file since 3772 was 3772, checked in by cbricaud, 11 years ago

add mutliplication by rn_pfac in qns computation for sf(jp_prec)%fnow(:,:,1) and sf(jp_snow)%fnow(:,:,1) ; see ticket 1048

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