New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
sbcblk_core.F90 in branches/2013/dev_LOCEAN_CMCC_INGV_MERC_UKMO_2013/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2013/dev_LOCEAN_CMCC_INGV_MERC_UKMO_2013/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 4245

Last change on this file since 4245 was 4245, checked in by cetlod, 11 years ago

dev_locean_cmcc_ingv_ukmo_merc : merge in the MERC_UKMO dev branch with trunk rev 4119

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