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

Last change on this file since 3432 was 3432, checked in by lguez, 5 years ago

Move the call to config_ocean_skin out of !$OMP MASTER in procedure
conf_phys. Definition of jcool, jwarm and rain_effect must be done by
all threads.

Use keywords in call to bulk_flux in procedure surf_ocean, for clarity.

Define CPP_KEY IN_LMDZ in makelmdz_fcm. This could be useful for any
external code used with LMDZ.

Move the Ocean_skin folder out of the LMDZ tree. This is more
convenient and clearer because Ocean_skin stays under Git control for
now. So we do not declare a phylmd/Ocean_skin folder in
"bld.cfg". Instead, we use the option -ext_src of makelmdz_fcm.

  • 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: 13.8 KB
Line 
1!
2! $Id: surf_ocean_mod.F90 3432 2019-01-18 17:29:40Z 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, t_int, s_int, ds_ns, dt_ns, dter, dser, tkt, tks, rf, &
23       taur)
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 sens_heat_rain_m, only: sens_heat_rain
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, DIMENSION(klon), INTENT(IN)        :: tsurf_in
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):: ds_ns(:) ! (knon)
79    ! "delta salinity near surface". Salinity variation in the
80    ! near-surface turbulent layer. That is subskin salinity minus
81    ! foundation salinity. In ppt.
82
83    REAL, intent(inout):: dt_ns(:) ! (knon)
84    ! "delta temperature near surface". Temperature variation in the
85    ! near-surface turbulent layer. That is subskin temperature
86    ! minus foundation temperature. (Can be negative.) In K.
87
88    ! Output variables
89    !******************************************************************************
90    REAL, DIMENSION(klon), INTENT(OUT)       :: z0m
91    !albedo SB >>>
92    !    REAL, DIMENSION(klon), INTENT(OUT)       :: alb1_new  ! new albedo in visible SW interval
93    !    REAL, DIMENSION(klon), INTENT(OUT)       :: alb2_new  ! new albedo in near IR interval
94    REAL, DIMENSION(6), INTENT(IN)          :: SFRWL
95    REAL, DIMENSION(klon,nsw), INTENT(OUT)       :: alb_dir_new,alb_dif_new
96    !albedo SB <<<     
97    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
98    REAL, INTENT(OUT):: tsurf_new(klon) ! sea surface temperature, in K
99    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l     
100    REAL, DIMENSION(klon), INTENT(OUT)       :: lmt_bils
101    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
102
103    REAL, intent(out):: t_int(:) ! (knon) interface temperature, in K
104    real, intent(out):: s_int(:) ! (knon) interface salinity, in ppt
105   
106    REAL, intent(out):: dter(:) ! (knon)
107    ! Temperature variation in the diffusive microlayer, that is
108    ! subskin temperature minus ocean-air interface temperature. In
109    ! K.
110
111    REAL, intent(out):: dser(:) ! (knon)
112    ! Temperature variation in the diffusive microlayer, that is
113    ! subskin temperature minus ocean-air interface temperature. In K.
114
115    REAL, intent(out):: tkt(:) ! (knon)
116    ! épaisseur (m) de la couche de diffusion thermique (microlayer)
117    ! cool skin thickness
118
119    REAL, intent(out):: tks(:) ! (knon)
120    ! épaisseur (m) de la couche de diffusion de masse (microlayer)
121
122    REAL, intent(out):: rf(:) ! (knon)
123    ! sensible heat flux at the surface due to rainfall, in  W m-2
124
125    REAL, intent(out):: taur(:) ! (knon)
126    ! momentum flux due to rain, in Pa
127
128    ! Local variables
129    !******************************************************************************
130    INTEGER               :: i, k
131    REAL                  :: tmp
132    REAL, PARAMETER       :: cepdu2=(0.1)**2
133    REAL, DIMENSION(klon) :: alb_eau, z0_lim
134    REAL, DIMENSION(klon) :: radsol
135    REAL, DIMENSION(klon) :: cdragq ! Cdrag pour l'evaporation
136    CHARACTER(len=20),PARAMETER :: modname="surf_ocean"
137    real rhoa(knon) ! density of moist air  (kg / m3)
138    real xlv(knon) ! chaleur latente d'évaporation (J / kg)
139    real precip_tot(knon) ! rain + snow
140    real S1(knon) ! salinity at depth_1, in ppt
141
142    ! End definition
143    !******************************************************************************
144
145
146    !******************************************************************************
147    ! Calculate total net radiance at surface
148    !
149    !******************************************************************************
150    radsol(1:klon) = 0.0 ! initialisation a priori inutile
151    radsol(1:knon) = swnet(1:knon) + lwnet(1:knon)
152
153    !******************************************************************************
154    ! Cdragq computed from cdrag
155    ! The difference comes only from a factor (f_z0qh_oce) on z0, so that
156    ! it can be computed inside surf_ocean
157    ! More complicated appraches may require the propagation through
158    ! pbl_surface of an independant cdragq variable.
159    !******************************************************************************
160
161    IF ( f_z0qh_oce .ne. 1.) THEN
162       ! Si on suit les formulations par exemple de Tessel, on
163       ! a z0h=0.4*nu/u*, z0q=0.62*nu/u*, d'ou f_z0qh_oce=0.62/0.4=1.55
164       cdragq(1:knon)=cdragh(1:knon)*                                      &
165            log(z1lay(1:knon)/z0h(1:knon))/log(z1lay(1:knon)/(f_z0qh_oce*z0h(1:knon)))
166    ELSE
167       cdragq(1:knon)=cdragh(1:knon)
168    ENDIF
169
170
171    !******************************************************************************
172    ! Switch according to type of ocean (couple, slab or forced)
173    !******************************************************************************
174    SELECT CASE(type_ocean)
175    CASE('couple')
176       CALL ocean_cpl_noice( &
177            swnet, lwnet, alb1, &
178            windsp, fder, &
179            itime, dtime, knon, knindex, &
180            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow,temp_air,spechum,&
181            AcoefH, AcoefQ, BcoefH, BcoefQ, &
182            AcoefU, AcoefV, BcoefU, BcoefV, &
183            ps, u1, v1, gustiness, &
184            radsol, snow, agesno, &
185            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
186            tsurf_new, dflux_s, dflux_l)
187
188    CASE('slab')
189       CALL ocean_slab_noice( &
190            itime, dtime, jour, knon, knindex, &
191            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, temp_air, spechum,&
192            AcoefH, AcoefQ, BcoefH, BcoefQ, &
193            AcoefU, AcoefV, BcoefU, BcoefV, &
194            ps, u1, v1, gustiness, tsurf_in, &
195            radsol, snow, &
196            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
197            tsurf_new, dflux_s, dflux_l, lmt_bils)
198
199    CASE('force')
200       CALL ocean_forced_noice( &
201            itime, dtime, jour, knon, knindex, &
202            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, &
203            temp_air, spechum, &
204            AcoefH, AcoefQ, BcoefH, BcoefQ, &
205            AcoefU, AcoefV, BcoefU, BcoefV, &
206            ps, u1, v1, gustiness, &
207            radsol, snow, agesno, &
208            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
209            tsurf_new, dflux_s, dflux_l)
210    END SELECT
211
212    !******************************************************************************
213    ! fcodron: compute lmt_bils  forced case (same as wfbils_oce / 1.-contfracatm)
214    !******************************************************************************
215    IF (type_ocean.NE.'slab') THEN
216       lmt_bils(1:klon)=0.
217       DO i=1,knon
218          lmt_bils(knindex(i))=(swnet(i)+lwnet(i)+fluxsens(i)+fluxlat(i)) &
219               *pctsrf(knindex(i),is_oce)/(1.-zmasq(knindex(i)))
220       END DO
221    END IF
222
223    !******************************************************************************
224    ! Calculate ocean surface albedo
225    !******************************************************************************
226    !albedo SB >>>
227    IF (iflag_albedo==0) THEN
228       !--old parametrizations of ocean surface albedo
229       !
230       IF (iflag_cycle_diurne.GE.1) THEN
231          !
232          CALL alboc_cd(rmu0,alb_eau)
233          !
234          !--ad-hoc correction for model radiative balance tuning
235          !--now outside alboc_cd routine
236          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
237          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.0),1.0)
238          !
239       ELSE
240          !
241          CALL alboc(REAL(jour),rlat,alb_eau)
242          !--ad-hoc correction for model radiative balance tuning
243          !--now outside alboc routine
244          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
245          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.04),0.60)
246          !
247       ENDIF
248       !
249       DO i =1, knon
250          DO  k=1,nsw
251             alb_dir_new(i,k) = alb_eau(knindex(i))
252          ENDDO
253       ENDDO
254       !IM 09122015 next line corresponds to the old way of doing in LMDZ5A/IPSLCM5A versions
255       !albedo for diffuse radiation is taken the same as for direct radiation
256       alb_dif_new(1:knon,:)=alb_dir_new(1:knon,:)
257       !IM 09122015 end
258       !
259    ELSE IF (iflag_albedo==1) THEN
260       !--new parametrization of ocean surface albedo by Sunghye Baek
261       !--albedo for direct and diffuse radiation are different
262       !
263       CALL ocean_albedo(knon,rmu0,knindex,windsp,SFRWL,alb_dir_new,alb_dif_new)
264       !
265       !--ad-hoc correction for model radiative balance tuning
266       alb_dir_new(1:knon,:) = fmagic*alb_dir_new(1:knon,:) + pmagic
267       alb_dif_new(1:knon,:) = fmagic*alb_dif_new(1:knon,:) + pmagic
268       alb_dir_new(1:knon,:)=MIN(MAX(alb_dir_new(1:knon,:),0.0),1.0)
269       alb_dif_new(1:knon,:)=MIN(MAX(alb_dif_new(1:knon,:),0.0),1.0)
270       !
271    ELSE IF (iflag_albedo==2) THEN
272       ! F. Codron albedo read from limit.nc
273       CALL limit_read_rug_alb(itime, dtime, jour,&
274            knon, knindex, z0_lim, alb_eau)
275       DO i =1, knon
276          DO  k=1,nsw
277             alb_dir_new(i,k) = alb_eau(i)
278          ENDDO
279       ENDDO
280       alb_dif_new=alb_dir_new
281    ENDIF
282    !albedo SB <<<
283
284    !******************************************************************************
285    ! Calculate the rugosity
286    !******************************************************************************
287    IF (iflag_z0_oce==0) THEN
288       DO i = 1, knon
289          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
290          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
291               +  0.11*14e-6 / SQRT(cdragm(i) * tmp)
292          z0m(i) = MAX(1.5e-05,z0m(i))
293       ENDDO
294       z0h(1:knon)=z0m(1:knon) ! En attendant mieux
295
296    ELSE IF (iflag_z0_oce==1) THEN
297       DO i = 1, knon
298          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
299          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
300               + 0.11*14e-6 / SQRT(cdragm(i) * tmp)
301          z0m(i) = MAX(1.5e-05,z0m(i))
302          z0h(i)=0.4*14e-6 / SQRT(cdragm(i) * tmp)
303       ENDDO
304    ELSE IF (iflag_z0_oce==-1) THEN
305       DO i = 1, knon
306          z0m(i) = z0min
307          z0h(i) = z0min
308       ENDDO
309    ELSE
310       CALL abort_physic(modname,'version non prevue',1)
311    ENDIF
312
313    rhoa = PS(:KNON) / (Rd * temp_air(:knon) * (1. + retv * spechum(:knon)))
314    xlv = rlvtt
315    precip_tot = precip_rain(:knon) + precip_snow(:knon)
316    rf =  sens_heat_rain(precip_tot, temp_air(:knon), spechum(:knon), rhoa, &
317         xlv, tsurf_new(:knon), ps(:knon))
318    s1 = 35.
319    call bulk_flux(tkt, tks, taur, dter, dser, t_int, s_int, ds_ns, dt_ns, &
320         u = windsp(:knon), t_ocean_1 = tsurf_new(:knon), s1 = s1, &
321         rain = precip_tot, hf = - fluxsens(:knon), hlb = - fluxlat(:knon), &
322         rnl = - lwnet(:knon), &
323         tau = sqrt(flux_u1(:knon)**2 + flux_v1(:knon)**2), rhoa = rhoa, &
324         xlv = xlv, rf = rf, dtime = dtime, rns = swnet(:knon))
325
326  END SUBROUTINE surf_ocean
327  !******************************************************************************
328  !
329END MODULE surf_ocean_mod
Note: See TracBrowser for help on using the repository browser.