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.
flxblk.F90 in trunk/NEMO/OPA_SRC/SBC – NEMO

source: trunk/NEMO/OPA_SRC/SBC/flxblk.F90 @ 84

Last change on this file since 84 was 84, checked in by opalod, 20 years ago

CT : UPDATE057 : # General syntax, alignement, comments corrections

# l_ctl alone replace the set (l_ctl .AND. lwp)
# Add of diagnostics which are activated when using l_ctl logical

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 42.1 KB
Line 
1MODULE flxblk
2   !!======================================================================
3   !!                       ***  MODULE  flxblk  ***
4   !! Ocean forcing:  bulk thermohaline forcing of the ocean (or ice)
5   !!=====================================================================
6#if defined key_flx_bulk_monthly || defined key_flx_bulk_daily
7   !!----------------------------------------------------------------------
8   !!   'key_flx_bulk_monthly'   or                            MONTHLY bulk
9   !!   'key_flx_bulk_daily'                                     DAILY bulk
10   !!----------------------------------------------------------------------
11   !!   flx_blk        : thermohaline fluxes from bulk
12   !!   flx_blk_albedo : albedo for ocean and ice (clear and overcast skies)
13   !!   flx_blk_declin : solar declinaison
14   !!----------------------------------------------------------------------
15   !! * Modules used
16   USE oce             ! ocean dynamics and tracers
17   USE dom_oce         ! ocean space and time domain
18   USE cpl_oce         ! ???
19   USE phycst          ! physical constants
20   USE daymod
21   USE blk_oce         ! bulk variables
22   USE flx_oce         ! forcings variables
23   USE ocfzpt          ! ???
24   USE in_out_manager
25   USE lbclnk
26
27   IMPLICIT NONE
28   PRIVATE
29
30   !! * Accessibility
31   PUBLIC flx_blk        ! routine called by flx.F90
32   PUBLIC flx_blk_albedo ! routine called by limflx.F90
33
34   !! * Module variables
35   INTEGER, PARAMETER  ::   &
36      jpintsr = 24          ! number of time step between sunrise and sunset
37      !                     ! uses for heat flux computation
38   LOGICAL ::   &
39      lbulk_init = .TRUE.   ! flag, bulk initialization done or not)
40
41   REAL(wp), DIMENSION(jpi,jpj) ::   &
42      stauc            ,  &   ! cloud optical depth
43      sbudyko   
44
45   !! * constants for bulk computation (flx_blk)
46   REAL(wp), DIMENSION(19)  ::  &
47      budyko                  ! BUDYKO's coefficient
48   ! BUDYKO's coefficient (cloudiness effect on LW radiation):
49   DATA budyko / 1.00, 0.98, 0.95, 0.92, 0.89, 0.86, 0.83, 0.80, 0.78, 0.75,  &
50      &          0.72, 0.69, 0.67, 0.64, 0.61, 0.58, 0.56, 0.53, 0.50 /
51   REAL(wp), DIMENSION(20)  :: &
52      tauco                  ! cloud optical depth coefficient
53   ! Cloud optical depth coefficient
54   DATA tauco / 6.6, 6.6, 7.0, 7.2, 7.1, 6.8, 6.5, 6.6, 7.1, 7.6,   &
55      &         6.6, 6.1, 5.6, 5.5, 5.8, 5.8, 5.6, 5.6, 5.6, 5.6 /
56   REAL(wp)  ::            &  ! constant values
57      zeps    = 1.e-20  ,  &
58      zeps0   = 1.e-13  ,  &
59      zeps1   = 1.e-06  ,  &
60      zzero   = 0.e0    ,  &
61      zone    = 1.0
62
63   !! * constants for albedo computation (flx_blk_albedo)
64   REAL(wp) ::   &
65      c1     = 0.05  ,     &   ! constants values
66      c2     = 0.10  ,     &
67      albice = 0.50  ,     &   !  albedo of melting ice in the arctic and antarctic (Shine & Hendersson-Sellers)
68      cgren  = 0.06  ,     &   !  correction of the snow or ice albedo to take into account
69                               !  effects of cloudiness (Grenfell & Perovich, 1984)
70      alphd  = 0.80  ,     &   !  coefficients for linear interpolation used to compute
71      alphdi = 0.72  ,     &   !  albedo between two extremes values (Pyane, 1972)
72      alphc  = 0.65  ,     &
73      zmue   = 0.40            !  cosine of local solar altitude
74
75   !! * constants for solar declinaison computation (flx_blk_declin)
76   REAL(wp) ::                &
77      a0  =  0.39507671   ,   &  ! coefficients
78      a1  = 22.85684301   ,   &
79      a2  = -0.38637317   ,   &
80      a3  =  0.15096535   ,   &
81      a4  = -0.00961411   ,   &
82      b1  = -4.29692073   ,   &
83      b2  =  0.05702074   ,   &
84      b3  = -0.09028607   ,   &
85      b4  =  0.00592797
86   !!----------------------------------------------------------------------
87   !!   OPA 9.0 , LODYC-IPSL  (2003)
88   !!----------------------------------------------------------------------
89
90CONTAINS
91
92   SUBROUTINE flx_blk( kt, psst )
93      !!---------------------------------------------------------------------------
94      !!                     ***  ROUTINE flx_blk  ***
95      !!                 
96      !!  ** Purpose :   Computation of the heat fluxes at ocean and snow/ice
97      !!       surface the solar heat at ocean and snow/ice surfaces and the
98      !!       sensitivity of total heat fluxes to the SST variations
99      !!         
100      !!  ** Method  :   The flux of heat at the ice and ocean surfaces are derived
101      !!       from semi-empirical ( or bulk ) formulae which relate the flux to
102      !!       the properties of the surface and of the lower atmosphere. Here, we
103      !!       follow the work of Oberhuber, 1988   
104      !!
105      !!  ** Action  :   call flx_blk_albedo to compute ocean and ice albedo
106      !!          computation of snow precipitation
107      !!          computation of solar flux at the ocean and ice surfaces
108      !!          computation of the long-wave radiation for the ocean and sea/ice
109      !!          computation of turbulent heat fluxes over water and ice
110      !!          computation of evaporation over water
111      !!          computation of total heat fluxes sensitivity over ice (dQ/dT)
112      !!          computation of latent heat flux sensitivity over ice (dQla/dT)
113      !!
114      !! History :
115      !!   8.0  !  97-06  (Louvain-La-Neuve)  Original code
116      !!   8.5  !  02-09  (C. Ethe , G. Madec )  F90: Free form and module
117      !!----------------------------------------------------------------------
118      !! * Arguments
119      INTEGER , INTENT( in  )  ::   kt        ! time step
120      REAL(wp), INTENT( in ), DIMENSION(jpi,jpj) ::   &
121         &                          psst      ! Sea Surface Temperature
122
123      !! * Local variables
124      INTEGER  ::             &
125         ji, jj, jt        ,  &  ! dummy loop indices
126         indaet            ,  &  !  = -1, 0, 1 for odd, normal and leap years resp.
127         iday              ,  &  ! integer part of day
128         indxb             ,  &  ! index for budyko coefficient
129         indxc                   ! index for cloud depth coefficient
130
131      REAL(wp)  ::            & 
132         zalat , zclat     ,  &  ! latitude in degrees
133         zmt1, zmt2, zmt3  ,  &  ! tempory air temperatures variables
134         ztatm3, ztatm4    ,  &  ! power 3 and 4 of air temperature
135         z4tatm3           ,  &  ! 4 * ztatm3
136         zcmue             ,  &  ! cosine of local solar altitude
137         zcmue2            ,  &  ! root of zcmue1
138         zscmue            ,  &  ! square-root of zcmue1
139         zpcmue            ,  &  ! zcmue1**1.4
140         zdecl             ,  &  ! solar declination
141         zsdecl , zcdecl   ,  &  ! sine and cosine of solar declination
142         zalbo             ,  &  ! albedo of sea-water
143         zalbi             ,  &  ! albedo of ice
144         ztamr             ,  &  ! air temperature minus triple point of water (rtt)
145         ztaevbk           ,  &  ! part of net longwave radiation
146         zevi , zevo       ,  &  ! vapour pressure of ice and ocean
147         zind1,zind2,zind3 ,  &  ! switch for testing the values of air temperature
148         zinda             ,  &  ! switch for testing the values of sea ice cover
149         zpis2             ,  &  ! pi / 2
150         z2pi                    ! 2 * pi
151
152      REAL(wp)  ::            & 
153         zxday             ,  &  ! day of year
154         ztsec             ,  &  ! time in seconds
155         zdist             ,  &  ! distance between the sun and the earth during the year
156         zdaycor           ,  &  ! corr. factor to take into account the variation of
157         !                       ! zday when calc. the solar rad.   
158         zesi, zeso        ,  &  ! vapour pressure of ice and ocean at saturation
159         zesi2             ,  &  ! root of zesi
160         zqsato            ,  &  ! humidity close to the ocean surface (at saturation)   
161         zqsati            ,  &  ! humidity close to the ice surface (at saturation)
162         zqsati2           ,  &  ! root of  zqsati
163         zdesidt           ,  &  ! derivative of zesi, function of ice temperature
164         zdteta            ,  &  ! diff. betw. sst and air temperature
165         zdeltaq           ,  &  ! diff. betw. spec. hum. and hum. close to the surface
166         ztvmoy, zobouks   ,  &  ! tempory scalars
167         zpsims, zpsihs, zpsils, zobouku, zxins, zpsimu ,  & 
168         zpsihu, zpsilu, zstab,zpsim, zpsih, zpsil      ,  & 
169         zvatmg, zcmn, zchn, zcln, zcmcmn, zdenum       ,  & 
170         zdtetar, ztvmoyr, zlxins, zcmn2, zchcm, zclcm , zcoef
171
172      REAL(wp)  ::            & 
173         zrhova            ,  &  ! air density per wind speed
174         zcsho , zcleo     ,  &  ! transfer coefficient over ocean
175         zcshi , zclei     ,  &  ! transfer coefficient over ice-free
176         zrhovacleo        ,  &  ! air density per wind speed per transfer coef.
177         zrhovacsho, zrhovaclei, zrhovacshi, & 
178         ztice3            ,  &  ! power 3 of ice temperature
179         zticemb, zticemb2 ,  &  ! tempory air temperatures variables
180         zdqlw_ice         ,  &  ! sensitivity of long-wave flux over ice
181         zdqsb_ice         ,  &  ! sensitivity of sensible heat flux over ice
182         zdqla_ice         ,  &  ! sensitivity of latent heat flux over ice
183         zdl, zdr                ! fractionnal part of latitude
184      REAL(wp), DIMENSION(jpi,jpj) :: & 
185         zpatm            ,  &   ! atmospheric pressure
186         zqatm            ,  &   ! specific humidity
187         zes              ,  &   ! vapour pressure at saturation
188         zev, zevsqr      ,  &   ! vapour pressure and his square-root
189         zrhoa            ,  &   ! air density
190         ztatm            ,  &   ! air temperature in Kelvins
191         zfrld            ,  &   ! fraction of sea ice cover
192         zcatm1           ,  &   ! fraction of cloud
193         zcldeff                 ! correction factor to account cloud effect
194      REAL(wp), DIMENSION(jpi,jpj) ::   & 
195         zalbocsd         ,  &   ! albedo of ocean
196         zalboos          ,  &   ! albedo of ocean under overcast sky
197         zalbics          ,  &   ! albedo of ice under clear sky
198         zalbios          ,  &   ! albedo of ice under overcast sky
199         zalbomu          ,  &   ! albedo of ocean when zcmue is 0.4
200         zqsro            ,  &   ! solar radiation over ocean
201         zqsrics          ,  &   ! solar radiation over ice under clear sky
202         zqsrios          ,  &   ! solar radiation over ice under overcast sky
203         zcldcor          ,  &   ! cloud correction
204         zlsrise, zlsset  ,  &   ! sunrise and sunset
205         zlmunoon         ,  &   ! local noon solar altitude
206         zdlha            ,  &   ! length of the ninstr segments of the solar day
207         zps              ,  &   ! sine of latitude per sine of solar decli.
208         zpc                     ! cosine of latitude per cosine of solar decli.
209
210      REAL(wp), DIMENSION(jpi,jpj) ::   & 
211         zqlw_oce         ,  &   ! long-wave heat flux over ocean
212         zqlw_ice         ,  &   ! long-wave heat flux over ice
213         zqla_oce         ,  &   ! latent heat flux over ocean
214         zqla_ice         ,  &   ! latent heat flux over ice
215         zqsb_oce         ,  &   ! sensible heat flux over ocean
216         zqsb_ice                ! sensible heat flux over ice
217 
218      REAL(wp), DIMENSION(jpi,jpj,jpintsr) ::    &
219         zlha             ,  &   ! local hour angle
220         zalbocs          ,  &   ! tempory var. of ocean albedo under clear sky
221         zsqsro           ,  &   ! tempory var. of solar rad. over ocean
222         zsqsrics         ,  &   ! temp. var. of solar rad. over ice under clear sky
223         zsqsrios                ! temp. var. of solar rad. over ice under overcast sky
224      !!---------------------------------------------------------------------
225
226      !---------------------
227      !   Initilization    !
228      !---------------------
229#if ! defined key_ice_lim
230      tn_ice(:,:) = psst(:,:)
231#endif
232
233      !  Determine cloud optical depths as a function of latitude (Chou et al., 1981).
234      !  and the correction factor for taking into account  the effect of clouds
235      !------------------------------------------------------
236      IF( lbulk_init ) THEN
237         DO jj = 1, jpj 
238            DO ji = 1 , jpi
239               zalat          = ( 90.e0 - ABS( gphit(ji,jj) ) ) /  5.e0
240               zclat          = ( 95.e0 -      gphit(ji,jj)   ) / 10.e0
241               indxb          = 1 + INT( zalat ) 
242               !  correction factor to account for the effect of clouds
243               sbudyko(ji,jj) = budyko(indxb) 
244               indxc          = 1 + INT( zclat ) 
245               zdl            = zclat - INT( zclat ) 
246               zdr            = 1.0 - zdl
247               stauc(ji,jj)   = zdr * tauco( indxc ) + zdl * tauco( indxc + 1 ) 
248            END DO
249         END DO
250         IF( nleapy == 1 ) THEN
251            yearday = 366.e0
252         ELSE IF( nleapy == 0 ) THEN
253            yearday = 365.e0
254         ELSEIF( nleapy == 30) THEN
255            yearday = 360.e0
256         ENDIF
257         lbulk_init = .FALSE.
258      ENDIF
259
260      zqlw_oce(:,:) = 0.e0
261      zqla_oce(:,:) = 0.e0
262      zqsb_oce(:,:) = 0.e0
263      zqlw_ice(:,:) = 0.e0
264      zqla_ice(:,:) = 0.e0
265      zqsb_ice(:,:) = 0.e0
266
267      zpis2       = rpi / 2.
268      z2pi        = 2. * rpi
269
270 !CDIR NOVERRCHK
271      DO jj = 1, jpj
272 !CDIR NOVERRCHK
273         DO ji = 1, jpi
274
275            ztatm (ji,jj) = 273.15 + tatm  (ji,jj)  !  air temperature in Kelvins
276            zcatm1(ji,jj) = 1.0    - catm  (ji,jj)  !  fractional cloud cover
277            zfrld (ji,jj) = 1.0    - freeze(ji,jj)  !  fractional sea ice cover
278            zpatm(ji,jj)  = 101000.               !  pressure
279     
280            !  Computation of air density, obtained from the equation of state for dry air.
281            zrhoa(ji,jj) = zpatm(ji,jj) / ( 287.04 * ztatm(ji,jj) )
282     
283            !  zes : Saturation water vapour
284            ztamr = ztatm(ji,jj) - rtt
285            zmt1  = SIGN( 17.269, ztamr )
286            zmt2  = SIGN( 21.875, ztamr )
287            zmt3  = SIGN( 28.200, -ztamr )
288            zes(ji,jj) = 611.0 * EXP (  ABS( ztamr ) * MIN ( zmt1, zmt2 )   &
289               &                      / ( ztatm(ji,jj) - 35.86  + MAX( zzero, zmt3 ) ) )
290
291            !  zev : vapour pressure  (hatm is relative humidity) 
292            zev(ji,jj)   = hatm(ji,jj) * zes(ji,jj) 
293            !  square-root of vapour pressure
294!CDIR NOVERRCHK
295            zevsqr(ji,jj) = SQRT( zev(ji,jj) * 0.01 )
296            !  zqapb  : specific humidity
297            zqatm(ji,jj) = 0.622 * zev(ji,jj) / ( zpatm(ji,jj) - 0.378 * zev(ji,jj) )
298
299
300            !----------------------------------------------------
301            !   Computation of snow precipitation (Ledley, 1985) |
302            !----------------------------------------------------
303
304            zmt1  =   253.0 - ztatm(ji,jj)
305            zmt2  = ( 272.0 - ztatm(ji,jj) ) / 38.0 
306            zmt3  = ( 281.0 - ztatm(ji,jj) ) / 18.0
307            zind1 = MAX( zzero, SIGN( zone, zmt1 ) )
308            zind2 = MAX( zzero, SIGN( zone, zmt2 ) )
309            zind3 = MAX( zzero, SIGN( zone, zmt3 ) )
310            ! total precipitation
311            tprecip(ji,jj) = watm(ji,jj)
312            ! solid  (snow) precipitation
313            sprecip(ji,jj) = tprecip(ji,jj) *       &
314               &             (           zind1      &
315               &               + ( 1.0 - zind1 ) * ( zind2 * ( 0.5 + zmt2 ) + ( 1.0 - zind2 ) *  zind3 * zmt3 ) ) 
316         END DO
317      END DO
318
319      !----------------------------------------------------------
320      !   Computation of albedo (need to calculates heat fluxes)|
321      !-----------------------------------------------------------
322     
323      CALL flx_blk_albedo( zalbios, zalboos, zalbics, zalbomu )
324
325      !-------------------------------------
326      !   Computation of solar irradiance. |
327      !----------------------------------------
328      indaet   = 1 
329      !  compution of the day of the year at which the fluxes have to be calculate
330      !--The date corresponds to the middle of the time step.
331      ztsec = ( 2 * INT ( ( kt - 1 ) / nfbulk ) + 1 ) * rdtbs2
332      zxday = MOD( ztsec , raass ) / rday + 1.0 
333
334      iday   = INT( zxday )
335      IF(l_ctl)   WRITE(numout,*) ' declin : iday ', iday, ' nfbulk= ', nfbulk
336      !   computation of the solar declination, his sine and his cosine
337      CALL flx_blk_declin( indaet, iday, zdecl )
338     
339      zdecl    = zdecl * rad
340      zsdecl   = SIN( zdecl )
341      zcdecl   = COS( zdecl )
342     
343      !  correction factor added for computation of shortwave flux to take into account the variation of
344      !  the distance between the sun and the earth during the year (Oberhuber 1988)
345      zdist    = zxday * z2pi / yearday
346      zdaycor  = 1.0 + 0.0013 * SIN( zdist ) + 0.0342 * COS( zdist )
347
348!CDIR NOVERRCHK
349      DO jj = 1, jpj
350!CDIR NOVERRCHK
351         DO ji = 1, jpi
352            !  product of sine of latitude and sine of solar declination
353            zps     (ji,jj) = SIN( gphit(ji,jj) * rad ) * zsdecl
354            !  product of cosine of latitude and cosine of solar declination
355            zpc     (ji,jj) = COS( gphit(ji,jj) * rad ) * zcdecl
356            !  computation of the both local time of sunrise and sunset
357            zlsrise (ji,jj) = ACOS( - SIGN( zone, zps(ji,jj) ) * MIN( zone, SIGN( zone, zps(ji,jj) )  &
358               &                     * ( zps(ji,jj) / zpc(ji,jj) ) ) ) 
359            zlsset  (ji,jj) = - zlsrise(ji,jj)
360            !  dividing the solar day into jpintsr segments of length zdlha
361            zdlha   (ji,jj) = ( zlsrise(ji,jj) - zlsset(ji,jj) ) / REAL( jpintsr )
362            !  computation of the local noon solar altitude
363            zlmunoon(ji,jj) = ASIN ( ( zps(ji,jj) + zpc(ji,jj) ) ) / rad
364           
365            !  cloud correction taken from Reed (1977) (imposed lower than 1)
366            zcldcor (ji,jj) = MIN( zone, ( 1.e0 - 0.62 * catm(ji,jj) + 0.0019 * zlmunoon(ji,jj) ) )
367         END DO
368      END DO
369
370         !  Computation of solar heat flux at each time of the day between sunrise and sunset.
371         !  We do this to a better optimisation of the code
372         !------------------------------------------------------       
373
374!CDIR NOVERRCHK   
375      DO jt = 1, jpintsr   
376         zcoef = FLOAT( jt ) - 0.5
377!CDIR NOVERRCHK     
378         DO jj = 1, jpj
379!CDIR NOVERRCHK
380            DO ji = 1, jpi
381               !  local hour angle
382               zlha (ji,jj,jt) = COS ( zlsrise(ji,jj) - zcoef * zdlha(ji,jj) )
383
384               ! cosine of local solar altitude
385               zcmue              = MAX ( zzero ,   zps(ji,jj) + zpc(ji,jj) * zlha (ji,jj,jt)  )
386               zcmue2             = 1368.0 * zcmue * zcmue
387               zscmue             = SQRT ( zcmue )
388               zpcmue             = zcmue**1.4
389               ! computation of sea-water albedo (Payne, 1972)
390               zalbocs(ji,jj,jt)  = 0.05 / ( 1.1 * zpcmue + 0.15 )
391               zalbo              = zcatm1(ji,jj) * zalbocs(ji,jj,jt) + catm(ji,jj) * zalboos(ji,jj)
392               ! solar heat flux absorbed at ocean surfaces (Zillman, 1972)
393               zevo               = zev(ji,jj) * 1.0e-05
394               zsqsro(ji,jj,jt)   =  ( 1.0 - zalbo ) * zdlha(ji,jj) * zcmue2                &
395                                   / ( ( zcmue + 2.7 ) * zevo + 1.085 * zcmue +  0.10 )
396               !  solar heat flux absorbed at sea/ice surfaces
397               !  Formulation of Shine and Crane, 1984 adapted for high albedo surfaces
398
399               !  For clear sky       
400               zevi               = zevo
401               zalbi              = zalbics(ji,jj)
402               zsqsrics(ji,jj,jt) =  ( 1.0 - zalbi ) * zdlha(ji,jj) * zcmue2                &
403                  &                / ( ( 1.0 + zcmue ) * zevi + 1.2 * zcmue + 0.0455 )
404
405               ! For overcast sky
406               zalbi              = zalbios(ji,jj)
407               zsqsrios(ji,jj,jt) = zdlha(ji,jj) *                                                           &
408                  &                 ( ( 53.5 + 1274.5 * zcmue      ) *  zscmue * ( 1.0 - 0.996  * zalbi ) )  &
409                  &                 / (  1.0 + 0.139  * stauc(ji,jj) *           ( 1.0 - 0.9435 * zalbi ) )
410            END DO
411         END DO
412      END DO
413
414
415      !  Computation of daily (between sunrise and sunset) solar heat flux absorbed
416      !  at the ocean and snow/ice surfaces.
417      !--------------------------------------------------------------------
418
419      zalbocsd(:,:) = 0.e0
420      zqsro   (:,:) = 0.e0
421      zqsrics (:,:) = 0.e0
422      zqsrios (:,:) = 0.e0
423
424      DO jt = 1, jpintsr 
425#   if defined key_vectopt_loop && ! defined key_autotasking
426         DO ji = 1, jpij 
427            zalbocsd(ji,1) = zalbocsd(ji,1) + zdlha   (ji,1) * zalbocs(ji,1,jt)   &
428               &                                             / MAX( 2.0 * zlsrise(ji,1) , zeps0 )
429            zqsro   (ji,1) = zqsro   (ji,1) + zsqsro  (ji,1,jt)
430            zqsrics (ji,1) = zqsrics (ji,1) + zsqsrics(ji,1,jt)
431            zqsrios (ji,1) = zqsrios (ji,1) + zsqsrios(ji,1,jt)
432         END DO
433#  else
434         DO jj = 1, jpj
435            DO ji = 1, jpi 
436               zalbocsd(ji,jj) = zalbocsd(ji,jj) + zdlha(ji,jj) * zalbocs(ji,jj,jt)   &
437                  &                                              / MAX( 2.0 * zlsrise(ji,jj) , zeps0 )
438               zqsro  (ji,jj)  = zqsro   (ji,jj) + zsqsro  (ji,jj,jt)
439               zqsrics(ji,jj)  = zqsrics (ji,jj) + zsqsrics(ji,jj,jt)
440               zqsrios(ji,jj)  = zqsrios (ji,jj) + zsqsrios(ji,jj,jt)
441            END DO
442         END DO
443#  endif
444      END DO
445
446      DO jj = 1, jpj
447         DO ji = 1, jpi 
448
449            !-------------------------------------------
450            !  Computation of shortwave radiation.
451            !-------------------------------------------
452
453            ! the solar heat flux absorbed at ocean and snow/ice surfaces
454            !------------------------------------------------------------
455
456            ! For ocean
457            qsr_oce(ji,jj) = srgamma * zcldcor(ji,jj) * zqsro(ji,jj) / z2pi
458            zinda          = SIGN( zone , -( -0.5 - zfrld(ji,jj) ) )
459            zinda          = 1.0 - MAX( zzero , zinda )
460            qsr_oce(ji,jj) = ( 1.- zinda ) * qsr_oce(ji,jj)
461
462            ! For snow/ice
463            qsr_ice(ji,jj) = ( zcatm1(ji,jj) * zqsrics(ji,jj) + catm(ji,jj) * zqsrios(ji,jj) ) / z2pi
464
465
466            ! Taking into account the ellipsity of the earth orbit
467            !-----------------------------------------------------
468
469            qsr_ice(ji,jj) = qsr_ice(ji,jj) * zdaycor
470            qsr_oce(ji,jj) = qsr_oce(ji,jj) * zdaycor
471
472            !  fraction of net shortwave radiation which is not absorbed in the
473            !  thin surface layer and penetrates inside the ice cover
474            !  ( Maykut and Untersteiner, 1971 ; Elbert anbd Curry, 1993 )
475            !------------------------------------------------------------------
476
477            fr1_i0(ji,jj) = 0.18  * zcatm1(ji,jj) + 0.35 * catm(ji,jj) 
478            fr2_i0(ji,jj) = 0.82  * zcatm1(ji,jj) + 0.65 * catm(ji,jj)
479
480            !---------------------------------------------------------------------------
481            !   Computation of long-wave radiation  ( Berliand 1952 ; all latitudes )
482            !---------------------------------------------------------------------------
483
484            ! tempory variables
485            ztatm3         = ztatm(ji,jj) * ztatm(ji,jj) * ztatm(ji,jj)
486            ztatm4         = ztatm3 * ztatm(ji,jj)
487            z4tatm3        = 4. * ztatm3
488            zcldeff(ji,jj) = 1.0 - sbudyko(ji,jj) * catm(ji,jj) * catm(ji,jj)   
489            ztaevbk        = ztatm4 * zcldeff(ji,jj) * ( 0.39 - 0.05 * zevsqr(ji,jj) ) 
490
491            !  Long-Wave for Ice
492            !----------------------
493            zqlw_ice(ji,jj) = - emic * stefan * ( ztaevbk + z4tatm3 * ( tn_ice(ji,jj) - ztatm(ji,jj) ) ) 
494
495            !  Long-Wave for Ocean
496            !-----------------------
497            zqlw_oce(ji,jj) = - emic * stefan * ( ztaevbk + z4tatm3 * ( psst  (ji,jj) - ztatm(ji,jj) ) ) 
498
499         END DO
500      END DO
501
502      !----------------------------------------
503      !  Computation of turbulent heat fluxes  ( Latent and sensible )
504      !----------------------------------------       
505      !CDIR NOVERRCHK
506      DO jj = 2 , jpjm1
507!ib   DO jj = 1 , jpj
508         !CDIR NOVERRCHK
509         DO  ji = 1, jpi
510
511            !  Turbulent heat fluxes over water
512            !----------------------------------
513
514            ! zeso     : vapour pressure at saturation of ocean
515            ! zqsato   : humidity close to the ocean surface (at saturation)
516            zeso          =  611.0 * EXP ( 17.2693884 * ( psst(ji,jj) - rtt ) * tmask(ji,jj,1) / ( psst(ji,jj) - 35.86 ) )
517            zqsato        = ( 0.622 * zeso ) / ( zpatm(ji,jj) - 0.378 * zeso )
518
519            !  Drag coefficients from Large and Pond (1981,1982)
520            !---------------------------------------------------
521   
522            !  Stability parameters
523            zdteta         = psst(ji,jj) - ztatm(ji,jj)
524            zdeltaq        = zqatm(ji,jj) - zqsato
525            ztvmoy         = ztatm(ji,jj) * ( 1. + 2.2e-3 * ztatm(ji,jj) * zqatm(ji,jj) )
526            zdenum         = MAX( vatm(ji,jj) * vatm(ji,jj) * ztvmoy, zeps )
527!i
528!i          if( zdenum == 0.e0 ) then
529!i               write(numout,*) 'flxblk  zdenum=0 ', ji,jj
530!i               zdenum = zeps
531!i          endif
532!i
533            zdtetar        = zdteta / zdenum
534            ztvmoyr        = ztvmoy * ztvmoy * zdeltaq / zdenum
535           
536            ! For stable atmospheric conditions
537            zobouks        = -70.0 * 10. * ( zdtetar + 3.2e-3 * ztvmoyr )
538            zobouks        = MAX( zzero , zobouks )
539            zpsims         = -7.0 * zobouks
540            zpsihs         =  zpsims
541            zpsils         =  zpsims
542
543            !  For unstable atmospheric conditions
544            zobouku        = -100.0 * 10.0 * ( zdtetar + 2.2e-3 * ztvmoyr )
545            zobouku        = MIN( zzero , zobouku )
546            zxins          = ( 1. - 16. * zobouku )**0.25
547            zlxins         = LOG( ( 1. + zxins * zxins ) / 2. )
548            zpsimu         = 2. * LOG( ( 1 + zxins ) / 2. )  + zlxins - 2. * ATAN( zxins ) + zpis2
549            zpsihu         = 2. * zlxins
550            zpsilu         = zpsihu
551
552            ! computation of intermediate values
553            zstab          = MAX( zzero , SIGN( zone , zdteta ) )
554            zpsim          = zstab * zpsimu + (1.0 - zstab ) * zpsims
555            zpsih          = zstab * zpsihu + (1.0 - zstab ) * zpsihs
556            zpsil          = zpsih
557           
558            zvatmg         = MAX( 0.032 * 1.5e-3 * vatm(ji,jj) * vatm(ji,jj) / grav, zeps )
559!i
560!!          if( zvatmg == 0.e0 ) then
561!!               write(numout,*) 'flxblk  zvatmg=0 ', ji,jj
562!!               zvatmg = zeps
563!!          endif
564!i
565
566            zcmn           = vkarmn / LOG ( 10. / zvatmg )
567            zcmn2          = zcmn * zcmn
568            zchn           = 0.0327 * zcmn
569            zcln           = 0.0346 * zcmn
570            zcmcmn         = 1 / ( 1 - zcmn * zpsim / vkarmn )
571            zchcm          = zcmcmn / ( 1 - zchn * zpsih / ( vkarmn * zcmn ) )
572            zclcm          = zchcm
573
574
575            !  Transfer cofficient zcsho and zcleo over ocean according to Large and Pond (1981,1982)
576            !--------------------------------------------------------------
577            zcsho          = zchn * zchcm
578            zcleo          = zcln * zclcm 
579
580
581            !   Computation of sensible and latent fluxes over Ocean
582            !----------------------------------------------------------------
583
584            !  computation of intermediate values
585            zrhova         = zrhoa(ji,jj) * vatm(ji,jj)
586            zrhovacsho     = zrhova * zcsho
587            zrhovacleo     = zrhova * zcleo
588
589            ! sensible heat flux
590            zqsb_oce(ji,jj) = zrhovacsho * 1004.0  * ( psst(ji,jj) - ztatm(ji,jj) ) 
591         
592            !  latent heat flux
593            zqla_oce(ji,jj) = zrhovacleo * 2.5e+06 * ( zqsato      - zqatm(ji,jj) )
594               
595            !  Calculate evaporation over water. (kg/m2/s)
596            !-------------------------------------------------
597            evap(ji,jj)    = zqla_oce(ji,jj) / cevap
598               
599               
600            !  Turbulent heat fluxes over snow/ice.
601            !--------------------------------------------------
602           
603            !  zesi     : vapour pressure at saturation of ice
604            !  zqsati   : humidity close to the ice surface (at saturation)
605            zesi           =  611.0 * EXP ( 21.8745587 * tmask(ji,jj,1)   &   ! tmask needed to avoid overflow in the exponential
606               &                                       * ( tn_ice(ji,jj) - rtt ) / ( tn_ice(ji,jj) - 7.66 ) )
607            zqsati         = ( 0.622 * zesi ) / ( zpatm(ji,jj) - 0.378 * zesi )
608               
609            !  computation of intermediate values
610            zticemb        = ( tn_ice(ji,jj) - 7.66 )
611            zticemb2       = zticemb * zticemb 
612            ztice3         = tn_ice(ji,jj) * tn_ice(ji,jj) * tn_ice(ji,jj)
613            zqsati2        = zqsati * zqsati
614            zesi2          = zesi * zesi
615            zdesidt        = zesi * ( 9.5 * LOG( 10.0 ) * ( rtt - 7.66 )  / zticemb2 )
616               
617            !  Transfer cofficient zcshi and zclei over ice. Assumed to be constant Parkinson 1979 ; Maykut 1982
618            !--------------------------------------------------------------------
619            zcshi          = 1.75e-03
620            zclei          = zcshi
621               
622            !  Computation of sensible and latent fluxes over ice
623            !----------------------------------------------------------------
624               
625            !  computation of intermediate values
626            zrhova          = zrhoa(ji,jj) * vatm(ji,jj)
627            zrhovacshi      = zrhova * zcshi * 2.834e+06
628            zrhovaclei      = zrhova * zclei * 1004.0
629           
630            !  sensible heat flux
631            zqsb_ice(ji,jj) = zrhovaclei * ( tn_ice(ji,jj) - ztatm(ji,jj) )
632           
633            !  latent heat flux
634            zqla_ice(ji,jj) = zrhovacshi * ( zqsati        - zqatm(ji,jj) )
635            qla_ice (ji,jj) = zqla_ice(ji,jj)
636             
637            !  Computation of sensitivity of non solar fluxes (dQ/dT)
638            !---------------------------------------------------------------
639               
640            !  computation of long-wave, sensible and latent flux sensitivity
641            zdqlw_ice       = 4.0 * emic * stefan * ztice3
642            zdqsb_ice       = zrhovaclei
643            zdqla_ice       = zrhovacshi * ( zdesidt * ( zqsati2 / zesi2 ) * ( zpatm(ji,jj) / 0.622 ) )         
644           
645            !  total non solar sensitivity
646            dqns_ice(ji,jj) = -( zdqlw_ice + zdqsb_ice + zdqla_ice ) 
647           
648            ! latent flux sensitivity
649            dqla_ice(ji,jj) = zdqla_ice
650           
651         END DO
652      END DO
653
654      ! total non solar heat flux over ice
655      qnsr_ice(:,:) = zqlw_ice(:,:) - zqsb_ice(:,:) - zqla_ice(:,:)
656      ! total non solar heat flux over water
657      qnsr_oce(:,:) = zqlw_oce(:,:) - zqsb_oce(:,:) - zqla_oce(:,:)
658
659      ! solid precipitations ( kg/m2/day -> kg/m2/s)
660      tprecip(:,:) = tprecip  (:,:) / rday 
661      ! snow  precipitations ( kg/m2/day -> kg/m2/s)
662      sprecip(:,:) = sprecip  (:,:) / rday 
663!i
664      CALL lbc_lnk( qsr_oce (:,:) , 'T', 1. )
665      CALL lbc_lnk( qnsr_oce(:,:) , 'T', 1. )
666      CALL lbc_lnk( qsr_ice (:,:) , 'T', 1. )
667      CALL lbc_lnk( qnsr_ice(:,:) , 'T', 1. )
668      CALL lbc_lnk( qla_ice (:,:) , 'T', 1. )
669      CALL lbc_lnk( dqns_ice(:,:) , 'T', 1. )
670      CALL lbc_lnk( dqla_ice(:,:) , 'T', 1. )
671      CALL lbc_lnk( fr1_i0  (:,:) , 'T', 1. )
672      CALL lbc_lnk( fr2_i0  (:,:) , 'T', 1. )
673      CALL lbc_lnk( tprecip (:,:) , 'T', 1. )
674      CALL lbc_lnk( sprecip (:,:) , 'T', 1. )
675      CALL lbc_lnk( evap    (:,:) , 'T', 1. )
676!i
677!i
678      qsr_oce (:,:) = qsr_oce (:,:)*tmask(:,:,1)
679      qnsr_oce(:,:) = qnsr_oce(:,:)*tmask(:,:,1)
680      qsr_ice (:,:) = qsr_ice (:,:)*tmask(:,:,1)
681      qnsr_ice(:,:) = qnsr_ice(:,:)*tmask(:,:,1)
682      qla_ice (:,:) = qla_ice (:,:)*tmask(:,:,1)
683      dqns_ice(:,:) = dqns_ice(:,:)*tmask(:,:,1)
684      dqla_ice(:,:) = dqla_ice(:,:)*tmask(:,:,1)
685      fr1_i0  (:,:) = fr1_i0  (:,:)*tmask(:,:,1)
686      fr2_i0  (:,:) = fr2_i0  (:,:)*tmask(:,:,1)
687      tprecip (:,:) = tprecip (:,:)*tmask(:,:,1)
688      sprecip (:,:) = sprecip (:,:)*tmask(:,:,1)
689      evap    (:,:) = evap    (:,:)*tmask(:,:,1)
690!i
691
692   END SUBROUTINE flx_blk
693
694
695#if defined key_ice_lim
696   !!----------------------------------------------------------------------
697   !!   'key_ice_lim'                                         LIM ice model
698   !!----------------------------------------------------------------------
699
700   SUBROUTINE flx_blk_albedo( palb , palcn , palbp , palcnp )
701      !!----------------------------------------------------------------------
702      !!               ***  ROUTINE flx_blk_albedo  ***
703      !!         
704      !! ** Purpose :   Computation of the albedo of the snow/ice system
705      !!      as well as the ocean one
706      !!       
707      !! ** Method  : - Computation of the albedo of snow or ice (choose the
708      !!      rignt one by a large number of tests
709      !!              - Computation of the albedo of the ocean
710      !!
711      !! References :
712      !!      Shine and Hendersson-Sellers 1985, JGR, 90(D1), 2243-2250.
713      !!
714      !! History :
715      !!  8.0   !  01-04  (LIM 1.0)
716      !!  8.5   !  03-07  (C. Ethe, G. Madec)  Optimization (old name:shine)
717      !!----------------------------------------------------------------------
718      !! * Modules used
719      USE ice                   ! ???
720
721      !! * Arguments
722      REAL(wp), DIMENSION(jpi,jpj), INTENT(out) ::  &
723         palb         ,     &    !  albedo of ice under overcast sky
724         palcn        ,     &    !  albedo of ocean under overcast sky
725         palbp        ,     &    !  albedo of ice under clear sky
726         palcnp                  !  albedo of ocean under clear sky
727
728      !! * Local variables
729      INTEGER ::    &
730         ji, jj                   ! dummy loop indices
731      REAL(wp) ::   & 
732         zmue14         ,     &   !  zmue**1.4
733         zalbpsnm       ,     &   !  albedo of ice under clear sky when snow is melting
734         zalbpsnf       ,     &   !  albedo of ice under clear sky when snow is freezing
735         zalbpsn        ,     &   !  albedo of snow/ice system when ice is coverd by snow
736         zalbpic        ,     &   !  albedo of snow/ice system when ice is free of snow
737         zithsn         ,     &   !  = 1 for hsn >= 0 ( ice is cov. by snow ) ; = 0 otherwise (ice is free of snow)
738         zitmlsn        ,     &   !  = 1 freezinz snow (sist >=rt0_snow) ; = 0 melting snow (sist<rt0_snow)
739         zihsc1         ,     &   !  = 1 hsn <= c1 ; = 0 hsn > c1
740         zihsc2                   !  = 1 hsn >= c2 ; = 0 hsn < c2
741      REAL(wp), DIMENSION(jpi,jpj) ::  &
742         zalbfz         ,     &   !  ( = alphdi for freezing ice ; = albice for melting ice )
743         zficeth                  !  function of ice thickness
744      LOGICAL , DIMENSION(jpi,jpj) ::  &
745         llmask
746      !!---------------------------------------------------------------------
747     
748      !-------------------------                                                             
749      !  Computation of  zficeth
750      !--------------------------
751     
752      llmask = (hsnif == 0.e0) .AND. ( sist >= rt0_ice )
753      WHERE ( llmask )   !  ice free of snow and melts
754         zalbfz = albice
755      ELSEWHERE                   
756         zalbfz = alphdi
757      END WHERE
758     
759      DO jj = 1, jpj
760         DO ji = 1, jpi
761            IF( hicif(ji,jj) > 1.5 ) THEN
762               zficeth(ji,jj) = zalbfz(ji,jj)
763            ELSEIF( hicif(ji,jj) > 1.0  .AND. hicif(ji,jj) <= 1.5 ) THEN
764               zficeth(ji,jj) = 0.472 + 2.0 * ( zalbfz(ji,jj) - 0.472 ) * ( hicif(ji,jj) - 1.0 )
765            ELSEIF( hicif(ji,jj) > 0.05 .AND. hicif(ji,jj) <= 1.0 ) THEN
766               zficeth(ji,jj) = 0.2467 + 0.7049 * hicif(ji,jj)                                &
767                  &                    - 0.8608 * hicif(ji,jj) * hicif(ji,jj)                 &
768                  &                    + 0.3812 * hicif(ji,jj) * hicif(ji,jj) * hicif (ji,jj)
769            ELSE
770               zficeth(ji,jj) = 0.1 + 3.6 * hicif(ji,jj) 
771            ENDIF
772         END DO
773      END DO
774     
775      !-----------------------------------------------
776      !    Computation of the snow/ice albedo system
777      !-------------------------- ---------------------
778     
779      !    Albedo of snow-ice for clear sky.
780      !-----------------------------------------------   
781      DO jj = 1, jpj
782         DO ji = 1, jpi
783            !  Case of ice covered by snow.             
784           
785            !  melting snow       
786            zihsc1       = 1.0 - MAX ( zzero , SIGN ( zone , - ( hsnif(ji,jj) - c1 ) ) )
787            zalbpsnm     = ( 1.0 - zihsc1 ) * ( zficeth(ji,jj) + hsnif(ji,jj) * ( alphd - zficeth(ji,jj) ) / c1 ) &
788               &                 + zihsc1   * alphd 
789            !  freezing snow               
790            zihsc2       = MAX ( zzero , SIGN ( zone , hsnif(ji,jj) - c2 ) )
791            zalbpsnf     = ( 1.0 - zihsc2 ) * ( albice + hsnif(ji,jj) * ( alphc - albice ) / c2 )                 &
792               &                 + zihsc2   * alphc 
793           
794            zitmlsn      =  MAX ( zzero , SIGN ( zone , sist(ji,jj) - rt0_snow ) )   
795            zalbpsn      =  zitmlsn * zalbpsnf + ( 1.0 - zitmlsn ) * zalbpsnm 
796           
797            !  Case of ice free of snow.
798            zalbpic      = zficeth(ji,jj) 
799           
800            ! albedo of the system   
801            zithsn       = 1.0 - MAX ( zzero , SIGN ( zone , - hsnif(ji,jj) ) )
802            palbp(ji,jj) =  zithsn * zalbpsn + ( 1.0 - zithsn ) *  zalbpic
803         END DO
804      END DO
805     
806      !    Albedo of snow-ice for overcast sky.
807      !---------------------------------------------- 
808      palb(:,:)   = palbp(:,:) + cgren                                           
809     
810      !--------------------------------------------
811      !    Computation of the albedo of the ocean
812      !-------------------------- -----------------                                                         
813     
814      !  Parameterization of Briegled and Ramanathan, 1982
815      zmue14      = zmue**1.4                                       
816      palcnp(:,:) = 0.05 / ( 1.1 * zmue14 + 0.15 )               
817     
818      !  Parameterization of Kondratyev, 1969 and Payne, 1972
819      palcn(:,:)  = 0.06                                                 
820     
821   END SUBROUTINE flx_blk_albedo
822
823# else
824   !!----------------------------------------------------------------------
825   !!   Default option :                                   NO sea-ice model
826   !!----------------------------------------------------------------------
827
828   SUBROUTINE flx_blk_albedo( palb , palcn , palbp , palcnp )
829      !!----------------------------------------------------------------------
830      !!               ***  ROUTINE flx_blk_albedo  ***
831      !!
832      !! ** Purpose :   Computation of the albedo of the snow/ice system
833      !!      as well as the ocean one
834      !!
835      !! ** Method  :   Computation of the albedo of snow or ice (choose the
836      !!      wright one by a large number of tests Computation of the albedo
837      !!      of the ocean
838      !!
839      !! History :
840      !!  8.0   !  01-04  (LIM 1.0)
841      !!  8.5   !  03-07  (C. Ethe, G. Madec)  Optimization (old name:shine)
842      !!----------------------------------------------------------------------
843      !! * Arguments
844      REAL(wp), DIMENSION(jpi,jpj), INTENT(out) ::  &
845         palb         ,     &    !  albedo of ice under overcast sky
846         palcn        ,     &    !  albedo of ocean under overcast sky
847         palbp        ,     &    !  albedo of ice under clear sky
848         palcnp                  !  albedo of ocean under clear sky
849
850      REAL(wp) ::   &
851         zmue14                 !  zmue**1.4
852      !!----------------------------------------------------------------------
853
854      !--------------------------------------------
855      !    Computation of the albedo of the ocean
856      !-------------------------- -----------------
857
858      !  Parameterization of Briegled and Ramanathan, 1982
859      zmue14      = zmue**1.4
860      palcnp(:,:) = 0.05 / ( 1.1 * zmue14 + 0.15 )
861
862      !  Parameterization of Kondratyev, 1969 and Payne, 1972
863      palcn(:,:)  = 0.06
864
865      palb (:,:)  = palcn(:,:)
866      palbp(:,:)  = palcnp(:,:)
867
868   END SUBROUTINE flx_blk_albedo
869
870#endif
871
872   SUBROUTINE flx_blk_declin( ky, kday, pdecl )
873      !!---------------------------------------------------------------------------
874      !!               ***  ROUTINE flx_blk_declin  ***
875      !!         
876      !! ** Purpose :   Computation of the solar declination for the day
877      !!         kday ( in decimal degrees ).
878      !!       
879      !! ** Method  :
880      !!
881      !! History :
882      !!         original    : 01-04 (LIM)
883      !!         addition    : 02-08 (C. Ethe, G. Madec)
884      !!---------------------------------------------------------------------
885      !! * Arguments
886      INTEGER, INTENT( in ) ::   &
887         ky  ,        &  ! = -1, 0, 1 for odd, normal and leap years resp.
888         kday            ! day of the year ( kday = 1 on january 1)
889      REAL(wp), INTENT(out) ::  &
890         pdecl            ! solar declination
891
892      !! * Local variables
893      REAL(wp) ::                & 
894         zday                ,   &  ! corresponding day of type year (cf. ky)
895         zp1, zp2, zp3, zp4         ! temporary scalars
896      !!---------------------------------------------------------------------
897     
898      zday = FLOAT( kday ) 
899     
900      IF( ky == 1 )  THEN
901         zday = zday - 0.5
902      ELSEIF( ky == 3 )  THEN
903         zday = zday - 1.
904      ELSE
905         zday = REAL( kday )
906      ENDIF
907     
908      zp1 = rpi * ( 2.0 * zday - 367.0 ) / yearday
909      zp2 = 2. * zp1
910      zp3 = 3. * zp1
911      zp4 = 4. * zp1
912     
913      pdecl  = a0                                                                      &
914         &   + a1 * COS( zp1 ) + a2 * COS( zp2 ) + a3 * COS( zp3 ) + a4 * COS( zp4 )   &
915         &   + b1 * SIN( zp1 ) + b2 * SIN( zp2 ) + b3 * SIN( zp3 ) + b4 * SIN( zp4 )
916     
917   END SUBROUTINE flx_blk_declin
918
919#else
920   !!----------------------------------------------------------------------
921   !!   Default option :           Empty module                     NO bulk
922   !!----------------------------------------------------------------------
923CONTAINS
924   SUBROUTINE flx_blk            ! Empty routine
925   END SUBROUTINE flx_blk
926#endif
927
928   !!======================================================================
929END MODULE flxblk
Note: See TracBrowser for help on using the repository browser.