source: LMDZ6/branches/Ocean_skin/libf/phylmd/surf_ocean_mod.F90 @ 3627

Last change on this file since 3627 was 3627, checked in by lguez, 4 years ago

If the ocean skin parameterization is working (passively or actively,
activate_ocean_skin >= 1) and we are coupled to the ocean then
receive bulk salinity of the surface layer of the ocean from the ocean
and feed it to procedure bulk_flux instead of the constant
value 35. If the ocean skin parameterization is working actively
(activate_ocean_skin == 2) and we are coupled to the ocean then send
ocean-air interface temperature to the ocean. We can only send
interface temperature from the previous time-step since communication
with the ocean is before the call to bulk_flux. In module cpl_mod,
define cpl_t_int with rank 1: no dimension for cpl_index because
t_int is only defined over ocean. New dummy argument sss of
procedures cpl_receive_ocean_fields and ocean_cpl_noice. New dummy
argument t_int of cpl_send_ocean_fields. In procedure
surf_ocean, rename local variable s1 to sss and give it the size
klon, which is required by the coupling machinery.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 14.3 KB
RevLine 
[781]1!
[2538]2! $Id: surf_ocean_mod.F90 3627 2020-02-03 13:27:46Z lguez $
3!
[781]4MODULE surf_ocean_mod
5
6  IMPLICIT NONE
7
8CONTAINS
[3395]9  !
10  !******************************************************************************
11  !
[888]12  SUBROUTINE surf_ocean(rlon, rlat, swnet, lwnet, alb1, &
[2243]13       windsp, rmu0, fder, tsurf_in, &
[781]14       itime, dtime, jour, knon, knindex, &
[2254]15       p1lay, z1lay, cdragh, cdragm, precip_rain, precip_snow, temp_air, spechum, &
[1067]16       AcoefH, AcoefQ, BcoefH, BcoefQ, &
17       AcoefU, AcoefV, BcoefU, BcoefV, &
[2240]18       ps, u1, v1, gustiness, rugoro, pctsrf, &
[888]19       snow, qsurf, agesno, &
[2243]20       z0m, z0h, SFRWL, alb_dir_new, alb_dif_new, evap, fluxsens, fluxlat, &
[1067]21       tsurf_new, dflux_s, dflux_l, lmt_bils, &
[3429]22       flux_u1, flux_v1, t_int, s_int, ds_ns, dt_ns, dter, dser, tkt, tks, rf, &
23       taur)
[1067]24
[3395]25    use albedo, only: alboc, alboc_cd
[3429]26    use bulk_flux_m, only: bulk_flux
[3395]27    USE dimphy, ONLY: klon, zmasq
28    USE surface_data, ONLY     : type_ocean
29    USE ocean_forced_mod, ONLY : ocean_forced_noice
30    USE ocean_slab_mod, ONLY   : ocean_slab_noice
31    USE ocean_cpl_mod, ONLY    : ocean_cpl_noice
32    USE indice_sol_mod, ONLY : nbsrf, is_oce
33    USE limit_read_mod
[3429]34    use sens_heat_rain_m, only: sens_heat_rain
[3458]35    use config_ocean_skin_m, only: activate_ocean_skin
[3395]36    !
37    ! This subroutine will make a call to ocean_XXX_noice according to the ocean mode (force,
38    ! slab or couple). The calculations of albedo and rugosity for the ocean surface are
39    ! done in here because they are identical for the different modes of ocean.
[1785]40
41
[793]42    INCLUDE "YOMCST.h"
[781]43
[2178]44    include "clesphys.h"
[2455]45    ! for cycle_diurne and for iflag_z0_oce==-1 (prescribed z0)
[2178]46
[3395]47    ! Input variables
48    !******************************************************************************
[781]49    INTEGER, INTENT(IN)                      :: itime, jour, knon
50    INTEGER, DIMENSION(klon), INTENT(IN)     :: knindex
51    REAL, INTENT(IN)                         :: dtime
52    REAL, DIMENSION(klon), INTENT(IN)        :: rlon, rlat
[888]53    REAL, DIMENSION(klon), INTENT(IN)        :: swnet  ! net shortwave radiation at surface 
54    REAL, DIMENSION(klon), INTENT(IN)        :: lwnet  ! net longwave radiation at surface 
55    REAL, DIMENSION(klon), INTENT(IN)        :: alb1   ! albedo in visible SW interval
[3429]56    REAL, DIMENSION(klon), INTENT(IN)        :: windsp ! wind at 10 m, in m s-1
[781]57    REAL, DIMENSION(klon), INTENT(IN)        :: rmu0 
58    REAL, DIMENSION(klon), INTENT(IN)        :: fder
[3627]59    REAL, INTENT(IN):: tsurf_in(klon) ! defined only for subscripts 1:knon
[2254]60    REAL, DIMENSION(klon), INTENT(IN)        :: p1lay,z1lay ! pression (Pa) et altitude (m) du premier niveau
[1067]61    REAL, DIMENSION(klon), INTENT(IN)        :: cdragh
62    REAL, DIMENSION(klon), INTENT(IN)        :: cdragm
[781]63    REAL, DIMENSION(klon), INTENT(IN)        :: precip_rain, precip_snow
64    REAL, DIMENSION(klon), INTENT(IN)        :: temp_air, spechum
[1067]65    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefH, AcoefQ, BcoefH, BcoefQ
66    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefU, AcoefV, BcoefU, BcoefV
[781]67    REAL, DIMENSION(klon), INTENT(IN)        :: ps
[2240]68    REAL, DIMENSION(klon), INTENT(IN)        :: u1, v1, gustiness
[781]69    REAL, DIMENSION(klon), INTENT(IN)        :: rugoro
70    REAL, DIMENSION(klon,nbsrf), INTENT(IN)  :: pctsrf
71
[3395]72    ! In/Output variables
73    !******************************************************************************
[888]74    REAL, DIMENSION(klon), INTENT(INOUT)     :: snow
75    REAL, DIMENSION(klon), INTENT(INOUT)     :: qsurf
[781]76    REAL, DIMENSION(klon), INTENT(INOUT)     :: agesno
[3389]77    REAL, DIMENSION(klon), INTENT(inOUT):: z0h
[781]78
[3429]79    REAL, intent(inout):: ds_ns(:) ! (knon)
80    ! "delta salinity near surface". Salinity variation in the
81    ! near-surface turbulent layer. That is subskin salinity minus
82    ! foundation salinity. In ppt.
83
84    REAL, intent(inout):: dt_ns(:) ! (knon)
85    ! "delta temperature near surface". Temperature variation in the
86    ! near-surface turbulent layer. That is subskin temperature
87    ! minus foundation temperature. (Can be negative.) In K.
88
[3395]89    ! Output variables
90    !******************************************************************************
[3389]91    REAL, DIMENSION(klon), INTENT(OUT)       :: z0m
[3395]92    !albedo SB >>>
93    !    REAL, DIMENSION(klon), INTENT(OUT)       :: alb1_new  ! new albedo in visible SW interval
94    !    REAL, DIMENSION(klon), INTENT(OUT)       :: alb2_new  ! new albedo in near IR interval
[2227]95    REAL, DIMENSION(6), INTENT(IN)          :: SFRWL
96    REAL, DIMENSION(klon,nsw), INTENT(OUT)       :: alb_dir_new,alb_dif_new
[3395]97    !albedo SB <<<     
[781]98    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
[3429]99    REAL, INTENT(OUT):: tsurf_new(klon) ! sea surface temperature, in K
[781]100    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l     
[996]101    REAL, DIMENSION(klon), INTENT(OUT)       :: lmt_bils
[1067]102    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
[781]103
[3627]104    REAL, intent(out):: t_int(:) ! (knon) ocean-air interface temperature, in K
105    real, intent(out):: s_int(:) ! (knon) ocean-air interface salinity, in ppt
[3458]106
[3429]107    REAL, intent(out):: dter(:) ! (knon)
108    ! Temperature variation in the diffusive microlayer, that is
[3601]109    ! ocean-air interface temperature minus subskin temperature. In
[3429]110    ! K.
111
112    REAL, intent(out):: dser(:) ! (knon)
[3627]113    ! Salinity variation in the diffusive microlayer, that is
114    ! ocean-air interface salinity minus subskin salinity. In ppt.
[3429]115
116    REAL, intent(out):: tkt(:) ! (knon)
117    ! épaisseur (m) de la couche de diffusion thermique (microlayer)
118    ! cool skin thickness
119
120    REAL, intent(out):: tks(:) ! (knon)
121    ! épaisseur (m) de la couche de diffusion de masse (microlayer)
122
123    REAL, intent(out):: rf(:) ! (knon)
124    ! sensible heat flux at the surface due to rainfall, in  W m-2
[3556]125    ! positive upward
[3429]126
127    REAL, intent(out):: taur(:) ! (knon)
128    ! momentum flux due to rain, in Pa
129
[3395]130    ! Local variables
131    !******************************************************************************
[2227]132    INTEGER               :: i, k
[1146]133    REAL                  :: tmp
134    REAL, PARAMETER       :: cepdu2=(0.1)**2
[3002]135    REAL, DIMENSION(klon) :: alb_eau, z0_lim
[888]136    REAL, DIMENSION(klon) :: radsol
[2254]137    REAL, DIMENSION(klon) :: cdragq ! Cdrag pour l'evaporation
[2391]138    CHARACTER(len=20),PARAMETER :: modname="surf_ocean"
[3429]139    real rhoa(knon) ! density of moist air  (kg / m3)
140    real xlv(knon) ! chaleur latente d'évaporation (J / kg)
141    real precip_tot(knon) ! rain + snow
[781]142
[3627]143    real sss(klon)
144    ! Bulk salinity of the surface layer of the ocean, in ppt. (Only
145    ! defined for subscripts 1:knon, but we have to declare it with
146    ! size klon because of the coupling machinery.)
147
[3395]148    ! End definition
149    !******************************************************************************
[888]150
151
[3395]152    !******************************************************************************
153    ! Calculate total net radiance at surface
154    !
155    !******************************************************************************
[2719]156    radsol(1:klon) = 0.0 ! initialisation a priori inutile
[888]157    radsol(1:knon) = swnet(1:knon) + lwnet(1:knon)
158
[3395]159    !******************************************************************************
160    ! Cdragq computed from cdrag
161    ! The difference comes only from a factor (f_z0qh_oce) on z0, so that
162    ! it can be computed inside surf_ocean
163    ! More complicated appraches may require the propagation through
164    ! pbl_surface of an independant cdragq variable.
165    !******************************************************************************
[2254]166
167    IF ( f_z0qh_oce .ne. 1.) THEN
[3395]168       ! Si on suit les formulations par exemple de Tessel, on
169       ! a z0h=0.4*nu/u*, z0q=0.62*nu/u*, d'ou f_z0qh_oce=0.62/0.4=1.55
[2718]170       cdragq(1:knon)=cdragh(1:knon)*                                      &
[3395]171            log(z1lay(1:knon)/z0h(1:knon))/log(z1lay(1:knon)/(f_z0qh_oce*z0h(1:knon)))
[2254]172    ELSE
[2718]173       cdragq(1:knon)=cdragh(1:knon)
[2254]174    ENDIF
175
[3002]176
[3395]177    !******************************************************************************
178    ! Switch according to type of ocean (couple, slab or forced)
179    !******************************************************************************
[996]180    SELECT CASE(type_ocean)
[781]181    CASE('couple')
[888]182       CALL ocean_cpl_noice( &
183            swnet, lwnet, alb1, &
[1067]184            windsp, fder, &
[781]185            itime, dtime, knon, knindex, &
[2254]186            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow,temp_air,spechum,&
[1067]187            AcoefH, AcoefQ, BcoefH, BcoefQ, &
188            AcoefU, AcoefV, BcoefU, BcoefV, &
[3463]189            ps, u1, v1, gustiness, tsurf_in, &
[888]190            radsol, snow, agesno, &
[1067]191            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
[3627]192            tsurf_new, dflux_s, dflux_l, sss)
[781]193
194    CASE('slab')
[888]195       CALL ocean_slab_noice( &
[996]196            itime, dtime, jour, knon, knindex, &
[2254]197            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, temp_air, spechum,&
[1067]198            AcoefH, AcoefQ, BcoefH, BcoefQ, &
199            AcoefU, AcoefV, BcoefU, BcoefV, &
[2240]200            ps, u1, v1, gustiness, tsurf_in, &
[2209]201            radsol, snow, &
[1067]202            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
[996]203            tsurf_new, dflux_s, dflux_l, lmt_bils)
[3395]204
[781]205    CASE('force')
[888]206       CALL ocean_forced_noice( &
207            itime, dtime, jour, knon, knindex, &
[2254]208            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, &
[781]209            temp_air, spechum, &
[1067]210            AcoefH, AcoefQ, BcoefH, BcoefQ, &
211            AcoefU, AcoefV, BcoefU, BcoefV, &
[3463]212            ps, u1, v1, gustiness, tsurf_in, &
[888]213            radsol, snow, agesno, &
[1067]214            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
[996]215            tsurf_new, dflux_s, dflux_l)
[781]216    END SELECT
217
[3395]218    !******************************************************************************
219    ! fcodron: compute lmt_bils  forced case (same as wfbils_oce / 1.-contfracatm)
220    !******************************************************************************
[2057]221    IF (type_ocean.NE.'slab') THEN
[3395]222       lmt_bils(1:klon)=0.
223       DO i=1,knon
224          lmt_bils(knindex(i))=(swnet(i)+lwnet(i)+fluxsens(i)+fluxlat(i)) &
225               *pctsrf(knindex(i),is_oce)/(1.-zmasq(knindex(i)))
226       END DO
[2057]227    END IF
228
[3395]229    !******************************************************************************
230    ! Calculate ocean surface albedo
231    !******************************************************************************
232    !albedo SB >>>
233    IF (iflag_albedo==0) THEN
234       !--old parametrizations of ocean surface albedo
235       !
236       IF (iflag_cycle_diurne.GE.1) THEN
237          !
238          CALL alboc_cd(rmu0,alb_eau)
239          !
240          !--ad-hoc correction for model radiative balance tuning
241          !--now outside alboc_cd routine
242          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
243          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.0),1.0)
244          !
245       ELSE
246          !
247          CALL alboc(REAL(jour),rlat,alb_eau)
248          !--ad-hoc correction for model radiative balance tuning
249          !--now outside alboc routine
250          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
251          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.04),0.60)
252          !
253       ENDIF
254       !
255       DO i =1, knon
256          DO  k=1,nsw
257             alb_dir_new(i,k) = alb_eau(knindex(i))
258          ENDDO
259       ENDDO
260       !IM 09122015 next line corresponds to the old way of doing in LMDZ5A/IPSLCM5A versions
261       !albedo for diffuse radiation is taken the same as for direct radiation
262       alb_dif_new(1:knon,:)=alb_dir_new(1:knon,:)
263       !IM 09122015 end
264       !
265    ELSE IF (iflag_albedo==1) THEN
266       !--new parametrization of ocean surface albedo by Sunghye Baek
267       !--albedo for direct and diffuse radiation are different
268       !
269       CALL ocean_albedo(knon,rmu0,knindex,windsp,SFRWL,alb_dir_new,alb_dif_new)
270       !
271       !--ad-hoc correction for model radiative balance tuning
272       alb_dir_new(1:knon,:) = fmagic*alb_dir_new(1:knon,:) + pmagic
273       alb_dif_new(1:knon,:) = fmagic*alb_dif_new(1:knon,:) + pmagic
274       alb_dir_new(1:knon,:)=MIN(MAX(alb_dir_new(1:knon,:),0.0),1.0)
275       alb_dif_new(1:knon,:)=MIN(MAX(alb_dif_new(1:knon,:),0.0),1.0)
276       !
277    ELSE IF (iflag_albedo==2) THEN
278       ! F. Codron albedo read from limit.nc
279       CALL limit_read_rug_alb(itime, dtime, jour,&
280            knon, knindex, z0_lim, alb_eau)
281       DO i =1, knon
282          DO  k=1,nsw
283             alb_dir_new(i,k) = alb_eau(i)
284          ENDDO
285       ENDDO
286       alb_dif_new=alb_dir_new
[781]287    ENDIF
[3395]288    !albedo SB <<<
[2227]289
[3395]290    !******************************************************************************
291    ! Calculate the rugosity
292    !******************************************************************************
293    IF (iflag_z0_oce==0) THEN
294       DO i = 1, knon
295          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
296          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
297               +  0.11*14e-6 / SQRT(cdragm(i) * tmp)
298          z0m(i) = MAX(1.5e-05,z0m(i))
299       ENDDO
300       z0h(1:knon)=z0m(1:knon) ! En attendant mieux
[2243]301
[3395]302    ELSE IF (iflag_z0_oce==1) THEN
303       DO i = 1, knon
304          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
305          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
306               + 0.11*14e-6 / SQRT(cdragm(i) * tmp)
307          z0m(i) = MAX(1.5e-05,z0m(i))
308          z0h(i)=0.4*14e-6 / SQRT(cdragm(i) * tmp)
309       ENDDO
310    ELSE IF (iflag_z0_oce==-1) THEN
311       DO i = 1, knon
312          z0m(i) = z0min
313          z0h(i) = z0min
314       ENDDO
315    ELSE
[2391]316       CALL abort_physic(modname,'version non prevue',1)
[3395]317    ENDIF
[3429]318
[3458]319    if (activate_ocean_skin >= 1) then
320       rhoa = PS(:KNON) / (Rd * temp_air(:knon) * (1. + retv * spechum(:knon)))
321       xlv = rlvtt
322       precip_tot = precip_rain(:knon) + precip_snow(:knon)
323       rf =  sens_heat_rain(precip_tot, temp_air(:knon), spechum(:knon), rhoa, &
[3463]324            xlv, tsurf_in(:knon), ps(:knon))
[3627]325       if (type_ocean /= 'couple') sss(:knon) = 35.
[3458]326       call bulk_flux(tkt, tks, taur, dter, dser, t_int, s_int, ds_ns, dt_ns, &
[3627]327            u = windsp(:knon), t_ocean_1 = tsurf_new(:knon), s1 = sss(:knon), &
[3458]328            rain = precip_tot, hf = - fluxsens(:knon), hlb = - fluxlat(:knon), &
329            rnl = - lwnet(:knon), &
330            tau = sqrt(flux_u1(:knon)**2 + flux_v1(:knon)**2), rhoa = rhoa, &
331            xlv = xlv, rf = rf, dtime = dtime, rns = swnet(:knon))
[3478]332       if (activate_ocean_skin == 2) tsurf_new(:knon) = t_int
[3458]333    end if
334   
[781]335  END SUBROUTINE surf_ocean
[3458]336  !****************************************************************************
[3395]337  !
[781]338END MODULE surf_ocean_mod
Note: See TracBrowser for help on using the repository browser.