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

source: branches/2012/dev_r3385_NOCS04_HAMF/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 3396

Last change on this file since 3396 was 3396, checked in by acc, 12 years ago

Branch: dev_r3385_NOCS04_HAMF; #665. Stage 1 of 2012 development: porting of changes on old development branch (2011/DEV_r1837_mass_heat_salt_fluxes) into new branch. Corrected a few errors on the way. This branch now compiles but is incomplete. Still missing LIM3 changes which must reside on a certain persons laptop somewhere

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