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

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

Store delta_sal instead of s_int

Revision analoguous to revision [3744], for salinity. Store as a state
variable and send to the ocean delta_sal, the difference between
ocean-air interface salinity and bulk salinity, instead of s_int,
the interface salinity.

So replace dummy argument s_int of procedure cpl_send_ocean_fields
by dummy argument delta_sal. Replace dummy argument s_int of
procedures ocean_cpl_noice and surf_ocean by dummy argument
delta_sal. Replace variable s_int of module phys_state_var_mod
by variable delta_sal. Rename local variable ys_int of procedure
pbl_surface to ydelta_sal. Set variable delta_sal of module
phys_state_var_mod to 0 for an appearing ocean fraction and a
missing startup field. Replace variable o_s_int of module
phys_output_ctrlout_mod by variable o_delta_sal.

Rename variables cpl_s_int and cpl_s_int_2D of module cpl_mod to
cpl_delta_sal and cpl_delta_sal_2D. Rename variable ids_s_int of
module oasis to ids_delta_sal. Change infosend(ids_delta_sal)%name
to "CODELSSS".

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