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

source: branches/2013/dev_r3853_CNRS9_ConfSetting/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 3875

Last change on this file since 3875 was 3875, checked in by clevy, 11 years ago

Configuration Setting/Step? 1, see ticket:#1074

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