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/2013/dev_LOCEAN_CMCC_INGV_2013/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2013/dev_LOCEAN_CMCC_INGV_2013/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 4230

Last change on this file since 4230 was 4230, checked in by cetlod, 10 years ago

dev_LOCEAN_CMCC_INGV_2013 : merge LOCEAN & CMCC_INGV branches, see ticket #1182

  • Property svn:keywords set to Id
File size: 54.4 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   USE lib_fortran     ! to use key_nosignedzero
45
46   IMPLICIT NONE
47   PRIVATE
48
49   PUBLIC   sbc_blk_core         ! routine called in sbcmod module
50   PUBLIC   blk_ice_core         ! routine called in sbc_ice_lim module
51   PUBLIC   turb_core_2z         ! routine calles in sbcblk_mfs module
52
53   INTEGER , PARAMETER ::   jpfld   = 9           ! maximum number of files to read
54   INTEGER , PARAMETER ::   jp_wndi = 1           ! index of 10m wind velocity (i-component) (m/s)    at T-point
55   INTEGER , PARAMETER ::   jp_wndj = 2           ! index of 10m wind velocity (j-component) (m/s)    at T-point
56   INTEGER , PARAMETER ::   jp_humi = 3           ! index of specific humidity               ( % )
57   INTEGER , PARAMETER ::   jp_qsr  = 4           ! index of solar heat                      (W/m2)
58   INTEGER , PARAMETER ::   jp_qlw  = 5           ! index of Long wave                       (W/m2)
59   INTEGER , PARAMETER ::   jp_tair = 6           ! index of 10m air temperature             (Kelvin)
60   INTEGER , PARAMETER ::   jp_prec = 7           ! index of total precipitation (rain+snow) (Kg/m2/s)
61   INTEGER , PARAMETER ::   jp_snow = 8           ! index of snow (solid prcipitation)       (kg/m2/s)
62   INTEGER , PARAMETER ::   jp_tdif = 9           ! index of tau diff associated to HF tau   (N/m2)   at T-point
63   
64   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf   ! structure of input fields (file informations, fields read)
65         
66   !                                             !!! CORE bulk parameters
67   REAL(wp), PARAMETER ::   rhoa =    1.22        ! air density
68   REAL(wp), PARAMETER ::   cpa  = 1000.5         ! specific heat of air
69   REAL(wp), PARAMETER ::   Lv   =    2.5e6       ! latent heat of vaporization
70   REAL(wp), PARAMETER ::   Ls   =    2.839e6     ! latent heat of sublimation
71   REAL(wp), PARAMETER ::   Stef =    5.67e-8     ! Stefan Boltzmann constant
72   REAL(wp), PARAMETER ::   Cice =    1.4e-3      ! iovi 1.63e-3     ! transfer coefficient over ice
73   REAL(wp), PARAMETER ::   albo =    0.066       ! ocean albedo assumed to be constant
74
75   !                                  !!* Namelist namsbc_core : CORE bulk parameters
76   LOGICAL  ::   ln_2m       ! logical flag for height of air temp. and hum
77   LOGICAL  ::   ln_taudif   ! logical flag to use the "mean of stress module - module of mean stress" data
78   REAL(wp) ::   rn_pfac     ! multiplication factor for precipitation
79   REAL(wp) ::   rn_efac     ! multiplication factor for evaporation (clem)
80   REAL(wp) ::   rn_vfac     ! multiplication factor for ice/ocean velocity in the calculation of wind stress (clem)
81
82   !! * Substitutions
83#  include "domzgr_substitute.h90"
84#  include "vectopt_loop_substitute.h90"
85   !!----------------------------------------------------------------------
86   !! NEMO/OPA 3.3 , NEMO-consortium (2010)
87   !! $Id$
88   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
89   !!----------------------------------------------------------------------
90CONTAINS
91
92   SUBROUTINE sbc_blk_core( kt )
93      !!---------------------------------------------------------------------
94      !!                    ***  ROUTINE sbc_blk_core  ***
95      !!                   
96      !! ** Purpose :   provide at each time step the surface ocean fluxes
97      !!      (momentum, heat, freshwater and runoff)
98      !!
99      !! ** Method  : (1) READ each fluxes in NetCDF files:
100      !!      the 10m wind velocity (i-component) (m/s)    at T-point
101      !!      the 10m wind velocity (j-component) (m/s)    at T-point
102      !!      the 10m or 2m specific humidity     ( % )
103      !!      the solar heat                      (W/m2)
104      !!      the Long wave                       (W/m2)
105      !!      the 10m or 2m air temperature       (Kelvin)
106      !!      the total precipitation (rain+snow) (Kg/m2/s)
107      !!      the snow (solid prcipitation)       (kg/m2/s)
108      !!      the tau diff associated to HF tau   (N/m2)   at T-point   (ln_taudif=T)
109      !!              (2) CALL blk_oce_core
110      !!
111      !!      C A U T I O N : never mask the surface stress fields
112      !!                      the stress is assumed to be in the (i,j) mesh referential
113      !!
114      !! ** Action  :   defined at each time-step at the air-sea interface
115      !!              - utau, vtau  i- and j-component of the wind stress
116      !!              - taum, wndm  wind stress and 10m wind modules at T-point
117      !!              - qns, qsr    non-solar and solar heat fluxes
118      !!              - emp         upward mass flux (evapo. - precip.)
119      !!              - sfx         salt flux due to freezing/melting (non-zero only if ice is present)
120      !!                            (set in limsbc(_2).F90)
121      !!----------------------------------------------------------------------
122      INTEGER, INTENT(in) ::   kt   ! ocean time step
123      !!
124      INTEGER  ::   ierror   ! return error code
125      INTEGER  ::   ifpr     ! dummy loop indice
126      INTEGER  ::   jfld     ! dummy loop arguments
127      INTEGER  ::   ios      ! Local integer output status for namelist read
128      !!
129      CHARACTER(len=100) ::  cn_dir   !   Root directory for location of core files
130      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i     ! array of namelist informations on the fields to read
131      TYPE(FLD_N) ::   sn_wndi, sn_wndj, sn_humi, sn_qsr       ! informations about the fields to be read
132      TYPE(FLD_N) ::   sn_qlw , sn_tair, sn_prec, sn_snow      !   "                                 "
133      TYPE(FLD_N) ::   sn_tdif                                 !   "                                 "
134      NAMELIST/namsbc_core/ cn_dir , ln_2m  , ln_taudif, rn_pfac, rn_efac, rn_vfac,  &
135         &                  sn_wndi, sn_wndj, sn_humi  , sn_qsr ,           &
136         &                  sn_qlw , sn_tair, sn_prec  , sn_snow, sn_tdif
137      !!---------------------------------------------------------------------
138
139      !                                         ! ====================== !
140      IF( kt == nit000 ) THEN                   !  First call kt=nit000  !
141         !                                      ! ====================== !
142         !
143         REWIND( numnam_ref )              ! Namelist namsbc_core in reference namelist : CORE bulk parameters
144         READ  ( numnam_ref, namsbc_core, IOSTAT = ios, ERR = 901)
145901      IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_core in reference namelist', lwp )
146
147         REWIND( numnam_cfg )              ! Namelist namsbc_core in configuration namelist : CORE bulk parameters
148         READ  ( numnam_cfg, namsbc_core, IOSTAT = ios, ERR = 902 )
149902      IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_core in configuration namelist', lwp )
150
151         WRITE ( numond, namsbc_core )
152         !                                         ! check: do we plan to use ln_dm2dc with non-daily forcing?
153         IF( ln_dm2dc .AND. sn_qsr%nfreqh /= 24 )   & 
154            &   CALL ctl_stop( 'sbc_blk_core: ln_dm2dc can be activated only with daily short-wave forcing' ) 
155         IF( ln_dm2dc .AND. sn_qsr%ln_tint ) THEN
156            CALL ctl_warn( 'sbc_blk_core: ln_dm2dc is taking care of the temporal interpolation of daily qsr',   &
157                 &         '              ==> We force time interpolation = .false. for qsr' )
158            sn_qsr%ln_tint = .false.
159         ENDIF
160         !                                         ! store namelist information in an array
161         slf_i(jp_wndi) = sn_wndi   ;   slf_i(jp_wndj) = sn_wndj
162         slf_i(jp_qsr ) = sn_qsr    ;   slf_i(jp_qlw ) = sn_qlw
163         slf_i(jp_tair) = sn_tair   ;   slf_i(jp_humi) = sn_humi
164         slf_i(jp_prec) = sn_prec   ;   slf_i(jp_snow) = sn_snow
165         slf_i(jp_tdif) = sn_tdif
166         !                 
167         lhftau = ln_taudif                        ! do we use HF tau information?
168         jfld = jpfld - COUNT( (/.NOT. lhftau/) )
169         !
170         ALLOCATE( sf(jfld), STAT=ierror )         ! set sf structure
171         IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_core: unable to allocate sf structure' )
172         DO ifpr= 1, jfld
173            ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) )
174            IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) )
175         END DO
176         !                                         ! fill sf with slf_i and control print
177         CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_core', 'flux formulation for ocean surface boundary condition', 'namsbc_core' )
178         !
179         sfx(:,:) = 0._wp                          ! salt flux; zero unless ice is present (computed in limsbc(_2).F90)
180         !
181      ENDIF
182
183      CALL fld_read( kt, nn_fsbc, sf )             ! input fields provided at the current time-step
184
185      !                                            ! compute the surface ocean fluxes using CORE bulk formulea
186      IF( MOD( kt - 1, nn_fsbc ) == 0 )   CALL blk_oce_core( kt, sf, sst_m, ssu_m, ssv_m )
187
188#if defined key_cice
189      IF( MOD( kt - 1, nn_fsbc ) == 0 )   THEN
190         qlw_ice(:,:,1)   = sf(jp_qlw)%fnow(:,:,1) 
191         qsr_ice(:,:,1)   = sf(jp_qsr)%fnow(:,:,1)
192         tatm_ice(:,:)    = sf(jp_tair)%fnow(:,:,1)         
193         qatm_ice(:,:)    = sf(jp_humi)%fnow(:,:,1)
194         tprecip(:,:)     = sf(jp_prec)%fnow(:,:,1) * rn_pfac
195         sprecip(:,:)     = sf(jp_snow)%fnow(:,:,1) * rn_pfac
196         wndi_ice(:,:)    = sf(jp_wndi)%fnow(:,:,1)
197         wndj_ice(:,:)    = sf(jp_wndj)%fnow(:,:,1)
198      ENDIF
199#endif
200      !
201   END SUBROUTINE sbc_blk_core
202   
203   
204   SUBROUTINE blk_oce_core( kt, sf, pst, pu, pv )
205      !!---------------------------------------------------------------------
206      !!                     ***  ROUTINE blk_core  ***
207      !!
208      !! ** Purpose :   provide the momentum, heat and freshwater fluxes at
209      !!      the ocean surface at each time step
210      !!
211      !! ** Method  :   CORE bulk formulea for the ocean using atmospheric
212      !!      fields read in sbc_read
213      !!
214      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
215      !!              - vtau    : j-component of the stress at V-point  (N/m2)
216      !!              - taum    : Wind stress module at T-point         (N/m2)
217      !!              - wndm    : Wind speed module at T-point          (m/s)
218      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
219      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
220      !!              - evap    : Evaporation over the ocean            (kg/m2/s)
221      !!              - emp     : evaporation minus precipitation       (kg/m2/s)
222      !!
223      !!  ** Nota  :   sf has to be a dummy argument for AGRIF on NEC
224      !!---------------------------------------------------------------------
225      INTEGER  , INTENT(in   )                 ::   kt    ! time step index
226      TYPE(fld), INTENT(inout), DIMENSION(:)   ::   sf    ! input data
227      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pst   ! surface temperature                      [Celcius]
228      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pu    ! surface current at U-point (i-component) [m/s]
229      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pv    ! surface current at V-point (j-component) [m/s]
230      !
231      INTEGER  ::   ji, jj               ! dummy loop indices
232      REAL(wp) ::   zcoef_qsatw, zztmp   ! local variable
233      REAL(wp), DIMENSION(:,:), POINTER ::   zwnd_i, zwnd_j    ! wind speed components at T-point
234      REAL(wp), DIMENSION(:,:), POINTER ::   zqsatw            ! specific humidity at pst
235      REAL(wp), DIMENSION(:,:), POINTER ::   zqlw, zqsb        ! long wave and sensible heat fluxes
236      REAL(wp), DIMENSION(:,:), POINTER ::   zqla, zevap       ! latent heat fluxes and evaporation
237      REAL(wp), DIMENSION(:,:), POINTER ::   Cd                ! transfer coefficient for momentum      (tau)
238      REAL(wp), DIMENSION(:,:), POINTER ::   Ch                ! transfer coefficient for sensible heat (Q_sens)
239      REAL(wp), DIMENSION(:,:), POINTER ::   Ce                ! tansfert coefficient for evaporation   (Q_lat)
240      REAL(wp), DIMENSION(:,:), POINTER ::   zst               ! surface temperature in Kelvin
241      REAL(wp), DIMENSION(:,:), POINTER ::   zt_zu             ! air temperature at wind speed height
242      REAL(wp), DIMENSION(:,:), POINTER ::   zq_zu             ! air spec. hum.  at wind speed height
243      !!---------------------------------------------------------------------
244      !
245      IF( nn_timing == 1 )  CALL timing_start('blk_oce_core')
246      !
247      CALL wrk_alloc( jpi,jpj, zwnd_i, zwnd_j, zqsatw, zqlw, zqsb, zqla, zevap )
248      CALL wrk_alloc( jpi,jpj, Cd, Ch, Ce, zst, zt_zu, zq_zu )
249      !
250      ! local scalars ( place there for vector optimisation purposes)
251      zcoef_qsatw = 0.98 * 640380. / rhoa
252     
253      zst(:,:) = pst(:,:) + rt0      ! convert SST from Celcius to Kelvin (and set minimum value far above 0 K)
254
255      ! ----------------------------------------------------------------------------- !
256      !      0   Wind components and module at T-point relative to the moving ocean   !
257      ! ----------------------------------------------------------------------------- !
258
259      ! ... components ( U10m - U_oce ) at T-point (unmasked)
260      zwnd_i(:,:) = 0.e0 
261      zwnd_j(:,:) = 0.e0
262#if defined key_cyclone
263# if defined key_vectopt_loop
264!CDIR COLLAPSE
265# endif
266      CALL wnd_cyc( kt, zwnd_i, zwnd_j )    ! add Manu !
267      DO jj = 2, jpjm1
268         DO ji = fs_2, fs_jpim1   ! vect. opt.
269            sf(jp_wndi)%fnow(ji,jj,1) = sf(jp_wndi)%fnow(ji,jj,1) + zwnd_i(ji,jj)
270            sf(jp_wndj)%fnow(ji,jj,1) = sf(jp_wndj)%fnow(ji,jj,1) + zwnd_j(ji,jj)
271         END DO
272      END DO
273#endif
274#if defined key_vectopt_loop
275!CDIR COLLAPSE
276#endif
277      DO jj = 2, jpjm1
278         DO ji = fs_2, fs_jpim1   ! vect. opt.
279            zwnd_i(ji,jj) = (  sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pu(ji-1,jj  ) + pu(ji,jj) )  )
280            zwnd_j(ji,jj) = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pv(ji  ,jj-1) + pv(ji,jj) )  )
281         END DO
282      END DO
283      CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. )
284      CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. )
285      ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked)
286!CDIR NOVERRCHK
287!CDIR COLLAPSE
288      wndm(:,:) = SQRT(  zwnd_i(:,:) * zwnd_i(:,:)   &
289         &             + zwnd_j(:,:) * zwnd_j(:,:)  ) * tmask(:,:,1)
290
291      ! ----------------------------------------------------------------------------- !
292      !      I   Radiative FLUXES                                                     !
293      ! ----------------------------------------------------------------------------- !
294   
295      ! ocean albedo assumed to be constant + modify now Qsr to include the diurnal cycle                    ! Short Wave
296      zztmp = 1. - albo
297      IF( ln_dm2dc ) THEN   ;   qsr(:,:) = zztmp * sbc_dcy( sf(jp_qsr)%fnow(:,:,1) ) * tmask(:,:,1)
298      ELSE                  ;   qsr(:,:) = zztmp *          sf(jp_qsr)%fnow(:,:,1)   * tmask(:,:,1)
299      ENDIF
300!CDIR COLLAPSE
301      zqlw(:,:) = (  sf(jp_qlw)%fnow(:,:,1) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:)  ) * tmask(:,:,1)   ! Long  Wave
302      ! ----------------------------------------------------------------------------- !
303      !     II    Turbulent FLUXES                                                    !
304      ! ----------------------------------------------------------------------------- !
305
306      ! ... specific humidity at SST and IST
307!CDIR NOVERRCHK
308!CDIR COLLAPSE
309      zqsatw(:,:) = zcoef_qsatw * EXP( -5107.4 / zst(:,:) ) 
310
311      ! ... NCAR Bulk formulae, computation of Cd, Ch, Ce at T-point :
312      IF( ln_2m ) THEN
313         !! If air temp. and spec. hum. are given at different height (2m) than wind (10m) :
314         CALL TURB_CORE_2Z(2.,10., zst   , sf(jp_tair)%fnow,         &
315            &                      zqsatw, sf(jp_humi)%fnow, wndm,   &
316            &                      Cd    , Ch              , Ce  ,   &
317            &                      zt_zu , zq_zu                   )
318      ELSE
319         !! If air temp. and spec. hum. are given at same height than wind (10m) :
320!gm bug?  at the compiling phase, add a copy in temporary arrays...  ==> check perf
321!         CALL TURB_CORE_1Z( 10., zst   (:,:), sf(jp_tair)%fnow(:,:),              &
322!            &                    zqsatw(:,:), sf(jp_humi)%fnow(:,:), wndm(:,:),   &
323!            &                    Cd    (:,:),             Ch  (:,:), Ce  (:,:)  )
324!gm bug
325! ARPDBG - this won't compile with gfortran. Fix but check performance
326! as per comment above.
327         CALL TURB_CORE_1Z( 10., zst   , sf(jp_tair)%fnow(:,:,1),       &
328            &                    zqsatw, sf(jp_humi)%fnow(:,:,1), wndm, &
329            &                    Cd    , Ch              , Ce    )
330      ENDIF
331
332      ! ... tau module, i and j component
333      DO jj = 1, jpj
334         DO ji = 1, jpi
335            zztmp = rhoa * wndm(ji,jj) * Cd(ji,jj)
336            taum  (ji,jj) = zztmp * wndm  (ji,jj)
337            zwnd_i(ji,jj) = zztmp * zwnd_i(ji,jj)
338            zwnd_j(ji,jj) = zztmp * zwnd_j(ji,jj)
339         END DO
340      END DO
341
342      ! ... add the HF tau contribution to the wind stress module?
343      IF( lhftau ) THEN 
344!CDIR COLLAPSE
345         taum(:,:) = taum(:,:) + sf(jp_tdif)%fnow(:,:,1)
346      ENDIF
347      CALL iom_put( "taum_oce", taum )   ! output wind stress module
348
349      ! ... utau, vtau at U- and V_points, resp.
350      !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines
351      DO jj = 1, jpjm1
352         DO ji = 1, fs_jpim1
353            utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj  ) )
354            vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji  ,jj+1) )
355         END DO
356      END DO
357      CALL lbc_lnk( utau(:,:), 'U', -1. )
358      CALL lbc_lnk( vtau(:,:), 'V', -1. )
359
360      !  Turbulent fluxes over ocean
361      ! -----------------------------
362      IF( ln_2m ) THEN
363         ! Values of temp. and hum. adjusted to 10m must be used instead of 2m values
364         zevap(:,:) = rn_efac * MAX( 0.e0, rhoa    *Ce(:,:)*( zqsatw(:,:) - zq_zu(:,:) ) * wndm(:,:) )   ! Evaporation
365         zqsb (:,:) =                      rhoa*cpa*Ch(:,:)*( zst   (:,:) - zt_zu(:,:) ) * wndm(:,:)     ! Sensible Heat
366      ELSE
367!CDIR COLLAPSE
368         zevap(:,:) = rn_efac * MAX( 0.e0, rhoa    *Ce(:,:)*( zqsatw(:,:) - sf(jp_humi)%fnow(:,:,1) ) * wndm(:,:) )   ! Evaporation
369!CDIR COLLAPSE
370         zqsb (:,:) =            rhoa*cpa*Ch(:,:)*( zst   (:,:) - sf(jp_tair)%fnow(:,:,1) ) * wndm(:,:)     ! Sensible Heat
371      ENDIF
372!CDIR COLLAPSE
373      zqla (:,:) = Lv * zevap(:,:)                                                              ! Latent Heat
374
375      IF(ln_ctl) THEN
376         CALL prt_ctl( tab2d_1=zqla  , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=Ce , clinfo2=' Ce  : ' )
377         CALL prt_ctl( tab2d_1=zqsb  , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=Ch , clinfo2=' Ch  : ' )
378         CALL prt_ctl( tab2d_1=zqlw  , clinfo1=' blk_oce_core: zqlw   : ', tab2d_2=qsr, clinfo2=' qsr : ' )
379         CALL prt_ctl( tab2d_1=zqsatw, clinfo1=' blk_oce_core: zqsatw : ', tab2d_2=zst, clinfo2=' zst : ' )
380         CALL prt_ctl( tab2d_1=utau  , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
381            &          tab2d_2=vtau  , clinfo2=              ' vtau : '  , mask2=vmask )
382         CALL prt_ctl( tab2d_1=wndm  , clinfo1=' blk_oce_core: wndm   : ')
383         CALL prt_ctl( tab2d_1=zst   , clinfo1=' blk_oce_core: zst    : ')
384      ENDIF
385       
386      ! ----------------------------------------------------------------------------- !
387      !     III    Total FLUXES                                                       !
388      ! ----------------------------------------------------------------------------- !
389     
390!CDIR COLLAPSE
391      emp (:,:) = (  zevap(:,:)                                          &   ! mass flux (evap. - precip.)
392         &         - sf(jp_prec)%fnow(:,:,1) * rn_pfac  ) * tmask(:,:,1)
393!CDIR COLLAPSE
394      qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)                                &   ! Downward Non Solar flux
395         &     - sf(jp_snow)%fnow(:,:,1) * rn_pfac * lfus                         &   ! remove latent melting heat for solid precip
396         &     - zevap(:,:) * pst(:,:) * rcp                                      &   ! remove evap heat content at SST
397         &     + ( sf(jp_prec)%fnow(:,:,1) - sf(jp_snow)%fnow(:,:,1) ) * rn_pfac  &   ! add liquid precip heat content at Tair
398         &     * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp                          &   
399         &     + sf(jp_snow)%fnow(:,:,1) * rn_pfac                                &   ! add solid  precip heat content at min(Tair,Tsnow)
400         &     * ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic 
401      !
402      CALL iom_put( "qlw_oce",   zqlw )                 ! output downward longwave heat over the ocean
403      CALL iom_put( "qsb_oce", - zqsb )                 ! output downward sensible heat over the ocean
404      CALL iom_put( "qla_oce", - zqla )                 ! output downward latent   heat over the ocean
405      CALL iom_put( "qhc_oce",   qns-zqlw+zqsb+zqla )   ! output downward heat content of E-P 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)  ) - rn_vfac * 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)  ) - rn_vfac * 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) - rn_vfac * 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) - rn_vfac * 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) - rn_vfac * 0.5 * ( pui(ji-1,jj  ) + pui(ji,jj) )  )
536               zwndj_t = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 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) ) - rn_vfac * 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) ) - rn_vfac * 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               ! iovino
575               IF( ff(ji,jj) .GT. 0._wp ) THEN
576                  z_qlw(ji,jj,jl) = ( 0.95 * sf(jp_qlw)%fnow(ji,jj,1) - Stef * pst(ji,jj,jl) * zst3 ) * tmask(ji,jj,1)
577               ELSE
578                  z_qlw(ji,jj,jl) = 0.95 * ( sf(jp_qlw)%fnow(ji,jj,1) - Stef * pst(ji,jj,jl) * zst3 ) * tmask(ji,jj,1)
579               ENDIF
580               ! lw sensitivity
581               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3                                               
582
583               ! ----------------------------!
584               !     II    Turbulent FLUXES  !
585               ! ----------------------------!
586
587               ! ... turbulent heat fluxes
588               ! Sensible Heat
589               z_qsb(ji,jj,jl) = rhoa * cpa * Cice * z_wnds_t(ji,jj) * ( pst(ji,jj,jl) - sf(jp_tair)%fnow(ji,jj,1) )
590               ! Latent Heat
591               p_qla(ji,jj,jl) = rn_efac * MAX( 0.e0, rhoa * Ls  * Cice * z_wnds_t(ji,jj)   &                           
592                  &                         * (  11637800. * EXP( -5897.8 / pst(ji,jj,jl) ) / rhoa - sf(jp_humi)%fnow(ji,jj,1)  ) )
593               ! Latent heat sensitivity for ice (Dqla/Dt)
594               p_dqla(ji,jj,jl) = rn_efac * zcoef_dqla * z_wnds_t(ji,jj) / ( zst2 ) * EXP( -5897.8 / pst(ji,jj,jl) )
595               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
596               z_dqsb(ji,jj,jl) = zcoef_dqsb * z_wnds_t(ji,jj)
597
598               ! ----------------------------!
599               !     III    Total FLUXES     !
600               ! ----------------------------!
601               ! Downward Non Solar flux
602               p_qns (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - p_qla (ji,jj,jl)     
603               ! Total non solar heat flux sensitivity for ice
604               p_dqns(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + p_dqla(ji,jj,jl) )   
605            END DO
606            !
607         END DO
608         !
609      END DO
610      !
611      !--------------------------------------------------------------------
612      ! FRACTIONs of net shortwave radiation which is not absorbed in the
613      ! thin surface layer and penetrates inside the ice cover
614      ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 )
615   
616!CDIR COLLAPSE
617      p_fr1(:,:) = ( 0.18 * ( 1.0 - zcoef_frca ) + 0.35 * zcoef_frca )
618!CDIR COLLAPSE
619      p_fr2(:,:) = ( 0.82 * ( 1.0 - zcoef_frca ) + 0.65 * zcoef_frca )
620       
621!CDIR COLLAPSE
622      p_tpr(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac      ! total precipitation [kg/m2/s]
623!CDIR COLLAPSE
624      p_spr(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac      ! solid precipitation [kg/m2/s]
625      CALL iom_put( 'snowpre', p_spr * 86400. )                  ! Snow precipitation
626      CALL iom_put( 'precip', p_tpr * 86400. )                   ! Total precipitation
627      !
628      IF(ln_ctl) THEN
629         CALL prt_ctl(tab3d_1=p_qla   , clinfo1=' blk_ice_core: p_qla  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb  : ', kdim=ijpl)
630         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice_core: z_qlw  : ', tab3d_2=p_dqla  , clinfo2=' p_dqla : ', kdim=ijpl)
631         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice_core: z_dqsb : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw : ', kdim=ijpl)
632         CALL prt_ctl(tab3d_1=p_dqns  , clinfo1=' blk_ice_core: p_dqns : ', tab3d_2=p_qsr   , clinfo2=' p_qsr  : ', kdim=ijpl)
633         CALL prt_ctl(tab3d_1=pst     , clinfo1=' blk_ice_core: pst    : ', tab3d_2=p_qns   , clinfo2=' p_qns  : ', kdim=ijpl)
634         CALL prt_ctl(tab2d_1=p_tpr   , clinfo1=' blk_ice_core: p_tpr  : ', tab2d_2=p_spr   , clinfo2=' p_spr  : ')
635         CALL prt_ctl(tab2d_1=p_taui  , clinfo1=' blk_ice_core: p_taui : ', tab2d_2=p_tauj  , clinfo2=' p_tauj : ')
636         CALL prt_ctl(tab2d_1=z_wnds_t, clinfo1=' blk_ice_core: z_wnds_t : ')
637      ENDIF
638
639      CALL wrk_dealloc( jpi,jpj, z_wnds_t )
640      CALL wrk_dealloc( jpi,jpj,pdim, z_qlw, z_qsb, z_dqlw, z_dqsb ) 
641      !
642      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_core')
643      !
644   END SUBROUTINE blk_ice_core
645 
646
647   SUBROUTINE TURB_CORE_1Z(zu, sst, T_a, q_sat, q_a,   &
648      &                        dU , Cd , Ch   , Ce   )
649      !!----------------------------------------------------------------------
650      !!                      ***  ROUTINE  turb_core  ***
651      !!
652      !! ** Purpose :   Computes turbulent transfert coefficients of surface
653      !!                fluxes according to Large & Yeager (2004)
654      !!
655      !! ** 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
656      !!      Momentum, Latent and sensible heat exchange coefficients
657      !!      Caution: this procedure should only be used in cases when air
658      !!      temperature (T_air), air specific humidity (q_air) and wind (dU)
659      !!      are provided at the same height 'zzu'!
660      !!
661      !! References :   Large & Yeager, 2004 : ???
662      !!----------------------------------------------------------------------
663      REAL(wp)                , INTENT(in   ) ::   zu      ! altitude of wind measurement       [m]
664      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   sst     ! sea surface temperature         [Kelvin]
665      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   T_a     ! potential air temperature       [Kelvin]
666      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   q_sat   ! sea surface specific humidity   [kg/kg]
667      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   q_a     ! specific air humidity           [kg/kg]
668      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   dU      ! wind module |U(zu)-U(0)|        [m/s]
669      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Cd      ! transfert coefficient for momentum       (tau)
670      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Ch      ! transfert coefficient for temperature (Q_sens)
671      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Ce      ! transfert coefficient for evaporation  (Q_lat)
672      !!
673      INTEGER :: j_itt
674      INTEGER , PARAMETER ::   nb_itt = 3
675      REAL(wp), PARAMETER ::   grav   = 9.8   ! gravity                       
676      REAL(wp), PARAMETER ::   kappa  = 0.4   ! von Karman s constant
677
678      REAL(wp), DIMENSION(:,:), POINTER  ::   dU10          ! dU                                   [m/s]
679      REAL(wp), DIMENSION(:,:), POINTER  ::   dT            ! air/sea temperature differeence      [K]
680      REAL(wp), DIMENSION(:,:), POINTER  ::   dq            ! air/sea humidity difference          [K]
681      REAL(wp), DIMENSION(:,:), POINTER  ::   Cd_n10        ! 10m neutral drag coefficient
682      REAL(wp), DIMENSION(:,:), POINTER  ::   Ce_n10        ! 10m neutral latent coefficient
683      REAL(wp), DIMENSION(:,:), POINTER  ::   Ch_n10        ! 10m neutral sensible coefficient
684      REAL(wp), DIMENSION(:,:), POINTER  ::   sqrt_Cd_n10   ! root square of Cd_n10
685      REAL(wp), DIMENSION(:,:), POINTER  ::   sqrt_Cd       ! root square of Cd
686      REAL(wp), DIMENSION(:,:), POINTER  ::   T_vpot        ! virtual potential temperature        [K]
687      REAL(wp), DIMENSION(:,:), POINTER  ::   T_star        ! turbulent scale of tem. fluct.
688      REAL(wp), DIMENSION(:,:), POINTER  ::   q_star        ! turbulent humidity of temp. fluct.
689      REAL(wp), DIMENSION(:,:), POINTER  ::   U_star        ! turb. scale of velocity fluct.
690      REAL(wp), DIMENSION(:,:), POINTER  ::   L             ! Monin-Obukov length                  [m]
691      REAL(wp), DIMENSION(:,:), POINTER  ::   zeta          ! stability parameter at height zu
692      REAL(wp), DIMENSION(:,:), POINTER  ::   U_n10         ! neutral wind velocity at 10m         [m]   
693      REAL(wp), DIMENSION(:,:), POINTER  ::   xlogt, xct, zpsi_h, zpsi_m
694     
695      INTEGER , DIMENSION(:,:), POINTER  ::   stab          ! 1st guess stability test integer
696      !!----------------------------------------------------------------------
697      !
698      IF( nn_timing == 1 )  CALL timing_start('TURB_CORE_1Z')
699      !
700      CALL wrk_alloc( jpi,jpj, stab )   ! integer
701      CALL wrk_alloc( jpi,jpj, dU10, dT, dq, Cd_n10, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd, L )
702      CALL wrk_alloc( jpi,jpj, T_vpot, T_star, q_star, U_star, zeta, U_n10, xlogt, xct, zpsi_h, zpsi_m )
703
704      !! * Start
705      !! Air/sea differences
706      dU10 = max(0.5, dU)     ! we don't want to fall under 0.5 m/s
707      dT = T_a - sst          ! assuming that T_a is allready the potential temp. at zzu
708      dq = q_a - q_sat
709      !!   
710      !! Virtual potential temperature
711      T_vpot = T_a*(1. + 0.608*q_a)
712      !!
713      !! Neutral Drag Coefficient
714      stab    = 0.5 + sign(0.5,dT)    ! stable : stab = 1 ; unstable : stab = 0
715      IF  ( ln_cdgw ) THEN
716        cdn_wave = cdn_wave - rsmall*(tmask(:,:,1)-1)
717        Cd_n10(:,:) =   cdn_wave
718      ELSE
719        Cd_n10  = 1E-3 * ( 2.7/dU10 + 0.142 + dU10/13.09 )    !   L & Y eq. (6a)
720      ENDIF
721      sqrt_Cd_n10 = sqrt(Cd_n10)
722      Ce_n10  = 1E-3 * ( 34.6 * sqrt_Cd_n10 )               !   L & Y eq. (6b)
723      Ch_n10  = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1-stab)) !   L & Y eq. (6c), (6d)
724      !!
725      !! Initializing transfert coefficients with their first guess neutral equivalents :
726      Cd = Cd_n10 ;  Ce = Ce_n10 ;  Ch = Ch_n10 ;  sqrt_Cd = sqrt(Cd)
727
728      !! * Now starting iteration loop
729      DO j_itt=1, nb_itt
730         !! Turbulent scales :
731         U_star  = sqrt_Cd*dU10                !   L & Y eq. (7a)
732         T_star  = Ch/sqrt_Cd*dT               !   L & Y eq. (7b)
733         q_star  = Ce/sqrt_Cd*dq               !   L & Y eq. (7c)
734
735         !! Estimate the Monin-Obukov length :
736         L  = (U_star**2)/( kappa*grav*(T_star/T_vpot + q_star/(q_a + 1./0.608)) )
737
738         !! Stability parameters :
739         zeta  = zu/L ;  zeta   = sign( min(abs(zeta),10.0), zeta )
740         zpsi_h  = psi_h(zeta)
741         zpsi_m  = psi_m(zeta)
742
743         IF  ( ln_cdgw ) THEN
744           sqrt_Cd=kappa/((kappa/sqrt_Cd_n10) - zpsi_m) ; Cd=sqrt_Cd*sqrt_Cd;
745         ELSE
746           !! Shifting the wind speed to 10m and neutral stability :
747           U_n10 = dU10*1./(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)) !  L & Y eq. (9a)
748
749           !! Updating the neutral 10m transfer coefficients :
750           Cd_n10  = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09)              !  L & Y eq. (6a)
751           sqrt_Cd_n10 = sqrt(Cd_n10)
752           Ce_n10  = 1E-3 * (34.6 * sqrt_Cd_n10)                           !  L & Y eq. (6b)
753           stab    = 0.5 + sign(0.5,zeta)
754           Ch_n10  = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab))           !  L & Y eq. (6c), (6d)
755
756           !! Shifting the neutral  10m transfer coefficients to ( zu , zeta ) :
757           !!
758           xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10) - zpsi_m)
759           Cd  = Cd_n10/(xct*xct) ;  sqrt_Cd = sqrt(Cd)
760         ENDIF
761         !!
762         xlogt = log(zu/10.) - zpsi_h
763         !!
764         xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10
765         Ch  = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct
766         !!
767         xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10
768         Ce  = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct
769         !!
770      END DO
771      !!
772      CALL wrk_dealloc( jpi,jpj, stab )   ! integer
773      CALL wrk_dealloc( jpi,jpj, dU10, dT, dq, Cd_n10, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd, L )
774      CALL wrk_dealloc( jpi,jpj, T_vpot, T_star, q_star, U_star, zeta, U_n10, xlogt, xct, zpsi_h, zpsi_m )
775      !
776      IF( nn_timing == 1 )  CALL timing_stop('TURB_CORE_1Z')
777      !
778    END SUBROUTINE TURB_CORE_1Z
779
780
781    SUBROUTINE TURB_CORE_2Z(zt, zu, sst, T_zt, q_sat, q_zt, dU, Cd, Ch, Ce, T_zu, q_zu)
782      !!----------------------------------------------------------------------
783      !!                      ***  ROUTINE  turb_core  ***
784      !!
785      !! ** Purpose :   Computes turbulent transfert coefficients of surface
786      !!                fluxes according to Large & Yeager (2004).
787      !!
788      !! ** 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
789      !!      Momentum, Latent and sensible heat exchange coefficients
790      !!      Caution: this procedure should only be used in cases when air
791      !!      temperature (T_air) and air specific humidity (q_air) are at 2m
792      !!      whereas wind (dU) is at 10m.
793      !!
794      !! References :   Large & Yeager, 2004 : ???
795      !!----------------------------------------------------------------------
796      REAL(wp), INTENT(in   )                     ::   zt       ! height for T_zt and q_zt                   [m]
797      REAL(wp), INTENT(in   )                     ::   zu       ! height for dU                              [m]
798      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   sst      ! sea surface temperature              [Kelvin]
799      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   T_zt     ! potential air temperature            [Kelvin]
800      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_sat    ! sea surface specific humidity         [kg/kg]
801      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_zt     ! specific air humidity                 [kg/kg]
802      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   dU       ! relative wind module |U(zu)-U(0)|       [m/s]
803      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cd       ! transfer coefficient for momentum         (tau)
804      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ch       ! transfer coefficient for sensible heat (Q_sens)
805      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ce       ! transfert coefficient for evaporation   (Q_lat)
806      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   T_zu     ! air temp. shifted at zu                     [K]
807      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   q_zu     ! spec. hum.  shifted at zu               [kg/kg]
808
809      INTEGER :: j_itt
810      INTEGER , PARAMETER :: nb_itt = 3              ! number of itterations
811      REAL(wp), PARAMETER ::   grav   = 9.8          ! gravity                       
812      REAL(wp), PARAMETER ::   kappa  = 0.4          ! von Karman's constant
813     
814      REAL(wp), DIMENSION(:,:), POINTER ::   dU10          ! dU                                [m/s]
815      REAL(wp), DIMENSION(:,:), POINTER ::   dT            ! air/sea temperature differeence   [K]
816      REAL(wp), DIMENSION(:,:), POINTER ::   dq            ! air/sea humidity difference       [K]
817      REAL(wp), DIMENSION(:,:), POINTER ::   Cd_n10        ! 10m neutral drag coefficient
818      REAL(wp), DIMENSION(:,:), POINTER ::   Ce_n10        ! 10m neutral latent coefficient
819      REAL(wp), DIMENSION(:,:), POINTER ::   Ch_n10        ! 10m neutral sensible coefficient
820      REAL(wp), DIMENSION(:,:), POINTER ::   sqrt_Cd_n10   ! root square of Cd_n10
821      REAL(wp), DIMENSION(:,:), POINTER ::   sqrt_Cd       ! root square of Cd
822      REAL(wp), DIMENSION(:,:), POINTER ::   T_vpot        ! virtual potential temperature        [K]
823      REAL(wp), DIMENSION(:,:), POINTER ::   T_star        ! turbulent scale of tem. fluct.
824      REAL(wp), DIMENSION(:,:), POINTER ::   q_star        ! turbulent humidity of temp. fluct.
825      REAL(wp), DIMENSION(:,:), POINTER ::   U_star        ! turb. scale of velocity fluct.
826      REAL(wp), DIMENSION(:,:), POINTER ::   L             ! Monin-Obukov length                  [m]
827      REAL(wp), DIMENSION(:,:), POINTER ::   zeta_u        ! stability parameter at height zu
828      REAL(wp), DIMENSION(:,:), POINTER ::   zeta_t        ! stability parameter at height zt
829      REAL(wp), DIMENSION(:,:), POINTER ::   U_n10         ! neutral wind velocity at 10m        [m]
830      REAL(wp), DIMENSION(:,:), POINTER ::   xlogt, xct, zpsi_hu, zpsi_ht, zpsi_m
831
832      INTEGER , DIMENSION(:,:), POINTER ::   stab          ! 1st stability test integer
833      !!----------------------------------------------------------------------
834      !
835      IF( nn_timing == 1 )  CALL timing_start('TURB_CORE_2Z')
836      !
837      CALL wrk_alloc( jpi,jpj, dU10, dT, dq, Cd_n10, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd, L )
838      CALL wrk_alloc( jpi,jpj, T_vpot, T_star, q_star, U_star, zeta_u, zeta_t, U_n10 )
839      CALL wrk_alloc( jpi,jpj, xlogt, xct, zpsi_hu, zpsi_ht, zpsi_m )
840      CALL wrk_alloc( jpi,jpj, stab )   ! interger
841
842      !! Initial air/sea differences
843      dU10 = max(0.5, dU)      !  we don't want to fall under 0.5 m/s
844      dT = T_zt - sst 
845      dq = q_zt - q_sat
846
847      !! Neutral Drag Coefficient :
848      stab = 0.5 + sign(0.5,dT)                 ! stab = 1  if dT > 0  -> STABLE
849      IF( ln_cdgw ) THEN
850        cdn_wave = cdn_wave - rsmall*(tmask(:,:,1)-1)
851        Cd_n10(:,:) =   cdn_wave
852      ELSE
853        Cd_n10  = 1E-3*( 2.7/dU10 + 0.142 + dU10/13.09 ) 
854      ENDIF
855      sqrt_Cd_n10 = sqrt(Cd_n10)
856      Ce_n10  = 1E-3*( 34.6 * sqrt_Cd_n10 )
857      Ch_n10  = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1 - stab))
858
859      !! Initializing transf. coeff. with their first guess neutral equivalents :
860      Cd = Cd_n10 ;  Ce = Ce_n10 ;  Ch = Ch_n10 ;  sqrt_Cd = sqrt(Cd)
861
862      !! Initializing z_u values with z_t values :
863      T_zu = T_zt ;  q_zu = q_zt
864
865      !!  * Now starting iteration loop
866      DO j_itt=1, nb_itt
867         dT = T_zu - sst ;  dq = q_zu - q_sat ! Updating air/sea differences
868         T_vpot = T_zu*(1. + 0.608*q_zu)    ! Updating virtual potential temperature at zu
869         U_star = sqrt_Cd*dU10                ! Updating turbulent scales :   (L & Y eq. (7))
870         T_star  = Ch/sqrt_Cd*dT              !
871         q_star  = Ce/sqrt_Cd*dq              !
872         !!
873         L = (U_star*U_star) &                ! Estimate the Monin-Obukov length at height zu
874              & / (kappa*grav/T_vpot*(T_star*(1.+0.608*q_zu) + 0.608*T_zu*q_star))
875         !! Stability parameters :
876         zeta_u  = zu/L  ;  zeta_u = sign( min(abs(zeta_u),10.0), zeta_u )
877         zeta_t  = zt/L  ;  zeta_t = sign( min(abs(zeta_t),10.0), zeta_t )
878         zpsi_hu = psi_h(zeta_u)
879         zpsi_ht = psi_h(zeta_t)
880         zpsi_m  = psi_m(zeta_u)
881         !!
882         !! Shifting the wind speed to 10m and neutral stability : (L & Y eq.(9a))
883!        U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u)))
884         U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m))
885         !!
886         !! Shifting temperature and humidity at zu :          (L & Y eq. (9b-9c))
887!        T_zu = T_zt - T_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t))
888         T_zu = T_zt - T_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht)
889!        q_zu = q_zt - q_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t))
890         q_zu = q_zt - q_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht)
891         !!
892         !! q_zu cannot have a negative value : forcing 0
893         stab = 0.5 + sign(0.5,q_zu) ;  q_zu = stab*q_zu
894         !!
895         IF( ln_cdgw ) THEN
896            sqrt_Cd=kappa/((kappa/sqrt_Cd_n10) - zpsi_m) ; Cd=sqrt_Cd*sqrt_Cd;
897         ELSE
898           !! Updating the neutral 10m transfer coefficients :
899           Cd_n10  = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09)    ! L & Y eq. (6a)
900           sqrt_Cd_n10 = sqrt(Cd_n10)
901           Ce_n10  = 1E-3 * (34.6 * sqrt_Cd_n10)                 ! L & Y eq. (6b)
902           stab    = 0.5 + sign(0.5,zeta_u)
903           Ch_n10  = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab)) ! L & Y eq. (6c-6d)
904           !!
905           !!
906           !! Shifting the neutral 10m transfer coefficients to (zu,zeta_u) :
907           xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)
908           Cd = Cd_n10/(xct*xct) ; sqrt_Cd = sqrt(Cd)
909         ENDIF
910         !!
911         xlogt = log(zu/10.) - zpsi_hu
912         !!
913         xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10
914         Ch  = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct
915         !!
916         xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10
917         Ce  = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct
918         !!
919         !!
920      END DO
921      !!
922      CALL wrk_dealloc( jpi,jpj, dU10, dT, dq, Cd_n10, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd, L )
923      CALL wrk_dealloc( jpi,jpj, T_vpot, T_star, q_star, U_star, zeta_u, zeta_t, U_n10 )
924      CALL wrk_dealloc( jpi,jpj, xlogt, xct, zpsi_hu, zpsi_ht, zpsi_m )
925      CALL wrk_dealloc( jpi,jpj, stab )   ! interger
926      !
927      IF( nn_timing == 1 )  CALL timing_stop('TURB_CORE_2Z')
928      !
929    END SUBROUTINE TURB_CORE_2Z
930
931
932    FUNCTION psi_m(zta)   !! Psis, L & Y eq. (8c), (8d), (8e)
933      !-------------------------------------------------------------------------------
934      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: zta
935
936      REAL(wp), PARAMETER :: pi = 3.141592653589793_wp
937      REAL(wp), DIMENSION(jpi,jpj)             :: psi_m
938      REAL(wp), DIMENSION(:,:), POINTER        :: X2, X, stabit
939      !-------------------------------------------------------------------------------
940
941      CALL wrk_alloc( jpi,jpj, X2, X, stabit )
942
943      X2 = sqrt(abs(1. - 16.*zta))  ;  X2 = max(X2 , 1.0) ;  X  = sqrt(X2)
944      stabit    = 0.5 + sign(0.5,zta)
945      psi_m = -5.*zta*stabit  &                                                          ! Stable
946         &    + (1. - stabit)*(2*log((1. + X)/2) + log((1. + X2)/2) - 2*atan(X) + pi/2)  ! Unstable
947
948      CALL wrk_dealloc( jpi,jpj, X2, X, stabit )
949      !
950    END FUNCTION psi_m
951
952
953    FUNCTION psi_h( zta )    !! Psis, L & Y eq. (8c), (8d), (8e)
954      !-------------------------------------------------------------------------------
955      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   zta
956      !
957      REAL(wp), DIMENSION(jpi,jpj)             ::   psi_h
958      REAL(wp), DIMENSION(:,:), POINTER        :: X2, X, stabit
959      !-------------------------------------------------------------------------------
960
961      CALL wrk_alloc( jpi,jpj, X2, X, stabit )
962
963      X2 = sqrt(abs(1. - 16.*zta))  ;  X2 = max(X2 , 1.) ;  X  = sqrt(X2)
964      stabit    = 0.5 + sign(0.5,zta)
965      psi_h = -5.*zta*stabit  &                                       ! Stable
966         &    + (1. - stabit)*(2.*log( (1. + X2)/2. ))                 ! Unstable
967
968      CALL wrk_dealloc( jpi,jpj, X2, X, stabit )
969      !
970    END FUNCTION psi_h
971 
972   !!======================================================================
973END MODULE sbcblk_core
Note: See TracBrowser for help on using the repository browser.