source: LMDZ6/trunk/libf/phylmd/surf_ocean_mod.F90 @ 5274

Last change on this file since 5274 was 5274, checked in by abarral, 9 hours ago

Replace yomcst.h by existing module

  • 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: 17.0 KB
Line 
1!
2! $Id: surf_ocean_mod.F90 5274 2024-10-25 13:41:23Z abarral $
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, precip_bs, 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, &
23       dt_ds, tkt, tks, taur, sss &
24#ifdef ISO
25        &       ,xtprecip_rain, xtprecip_snow,xtspechum,Roce, &
26        &       xtsnow,xtevap,h1 &
27#endif               
28        &       )
29
30    use albedo, only: alboc, alboc_cd
31    use bulk_flux_m, only: bulk_flux
32    USE dimphy, ONLY: klon, zmasq
33    USE surface_data, ONLY     : type_ocean
34    USE ocean_forced_mod, ONLY : ocean_forced_noice
35    USE ocean_slab_mod, ONLY   : ocean_slab_noice
36    USE ocean_cpl_mod, ONLY    : ocean_cpl_noice
37    USE indice_sol_mod, ONLY : nbsrf, is_oce
38#ifdef ISO
39    USE infotrac_phy, ONLY : ntraciso=>ntiso,niso
40#ifdef ISOVERIF
41    USE isotopes_mod, ONLY: iso_eau,ridicule
42    USE isotopes_verif_mod
43#endif
44#endif
45    USE yomcst_mod_h, ONLY: RPI, RCLUM, RHPLA, RKBOL, RNAVO                   &
46          , RDAY, REA, REPSM, RSIYEA, RSIDAY, ROMEGA                  &
47          , R_ecc, R_peri, R_incl                                      &
48          , RA, RG, R1SA                                         &
49          , RSIGMA                                                     &
50          , R, RMD, RMV, RD, RV, RCPD                    &
51          , RMO3, RMCO2, RMC, RMCH4, RMN2O, RMCFC11, RMCFC12        &
52          , RCPV, RCVD, RCVV, RKAPPA, RETV, eps_w                    &
53          , RCW, RCS                                                 &
54          , RLVTT, RLSTT, RLMLT, RTT, RATM                           &
55          , RESTT, RALPW, RBETW, RGAMW, RALPS, RBETS, RGAMS            &
56          , RALPD, RBETD, RGAMD
57USE limit_read_mod
58    USE config_ocean_skin_m, ONLY: activate_ocean_skin
59    !
60    ! This subroutine will make a call to ocean_XXX_noice according to the ocean mode (force,
61    ! slab or couple). The calculations of albedo and rugosity for the ocean surface are
62    ! done in here because they are identical for the different modes of ocean.
63
64
65
66
67    include "clesphys.h"
68    ! for cycle_diurne and for iflag_z0_oce==-1 (prescribed z0)
69
70    ! Input variables
71    !******************************************************************************
72    INTEGER, INTENT(IN)                      :: itime, jour, knon
73    INTEGER, DIMENSION(klon), INTENT(IN)     :: knindex
74    REAL, INTENT(IN)                         :: dtime
75    REAL, DIMENSION(klon), INTENT(IN)        :: rlon, rlat
76    REAL, DIMENSION(klon), INTENT(IN)        :: swnet  ! net shortwave radiation at surface 
77    REAL, DIMENSION(klon), INTENT(IN)        :: lwnet  ! net longwave radiation at surface 
78    REAL, DIMENSION(klon), INTENT(IN)        :: alb1   ! albedo in visible SW interval
79    REAL, DIMENSION(klon), INTENT(IN)        :: windsp ! wind at 10 m, in m s-1
80    REAL, DIMENSION(klon), INTENT(IN)        :: rmu0 
81    REAL, DIMENSION(klon), INTENT(IN)        :: fder
82    REAL, DIMENSION(klon), INTENT(IN)        :: tsurf_in     ! defined only for subscripts 1:knon
83    REAL, DIMENSION(klon), INTENT(IN)        :: p1lay,z1lay ! pression (Pa) et altitude (m) du premier niveau
84    REAL, DIMENSION(klon), INTENT(IN)        :: cdragh
85    REAL, DIMENSION(klon), INTENT(IN)        :: cdragm
86    REAL, DIMENSION(klon), INTENT(IN)        :: precip_rain, precip_snow, precip_bs
87    REAL, DIMENSION(klon), INTENT(IN)        :: temp_air, spechum
88    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefH, AcoefQ, BcoefH, BcoefQ
89    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefU, AcoefV, BcoefU, BcoefV
90    REAL, DIMENSION(klon), INTENT(IN)        :: ps
91    REAL, DIMENSION(klon), INTENT(IN)        :: u1, v1, gustiness
92    REAL, DIMENSION(klon), INTENT(IN)        :: rugoro
93    REAL, DIMENSION(klon,nbsrf), INTENT(IN)  :: pctsrf
94#ifdef ISO
95    REAL, DIMENSION(ntraciso,klon), INTENT(IN) :: xtprecip_rain, xtprecip_snow
96    REAL, DIMENSION(ntraciso,klon), INTENT(IN) :: xtspechum
97#endif
98
99    ! In/Output variables
100    !******************************************************************************
101    REAL, DIMENSION(klon), INTENT(INOUT)     :: snow
102    REAL, DIMENSION(klon), INTENT(INOUT)     :: qsurf
103    REAL, DIMENSION(klon), INTENT(INOUT)     :: agesno
104    REAL, DIMENSION(klon), INTENT(inOUT)     :: z0h
105#ifdef ISO
106    REAL, DIMENSION(niso,klon), INTENT(IN)   :: xtsnow
107    REAL, DIMENSION(niso,klon), INTENT(INOUT):: Roce 
108#endif
109
110    REAL, intent(inout):: delta_sst(:) ! (knon)
111    ! Ocean-air interface temperature minus bulk SST, in K. Defined
112    ! only if activate_ocean_skin >= 1.
113
114    real, intent(inout):: delta_sal(:) ! (knon)
115    ! Ocean-air interface salinity minus bulk salinity, in ppt. Defined
116    ! only if activate_ocean_skin >= 1.
117
118    REAL, intent(inout):: ds_ns(:) ! (knon)
119    ! "delta salinity near surface". Salinity variation in the
120    ! near-surface turbulent layer. That is subskin salinity minus
121    ! foundation salinity. In ppt.
122
123    REAL, intent(inout):: dt_ns(:) ! (knon)
124    ! "delta temperature near surface". Temperature variation in the
125    ! near-surface turbulent layer. That is subskin temperature
126    ! minus foundation temperature. (Can be negative.) In K.
127
128    REAL, intent(inout):: dter(:) ! (knon)
129    ! Temperature variation in the diffusive microlayer, that is
130    ! ocean-air interface temperature minus subskin temperature. In
131    ! K.
132
133    REAL, intent(inout):: dser(:) ! (knon)
134    ! Salinity variation in the diffusive microlayer, that is
135    ! ocean-air interface salinity minus subskin salinity. In ppt.
136
137    real, intent(inout):: dt_ds(:) ! (knon)
138    ! (tks / tkt) * dTer, in K
139
140    ! Output variables
141    !**************************************************************************
142    REAL, DIMENSION(klon), INTENT(OUT)       :: z0m
143    !albedo SB >>>
144    !    REAL, DIMENSION(klon), INTENT(OUT)  :: alb1_new  ! new albedo in visible SW interval
145    !    REAL, DIMENSION(klon), INTENT(OUT)  :: alb2_new  ! new albedo in near IR interval
146    REAL, DIMENSION(6), INTENT(IN)           :: SFRWL
147    REAL, DIMENSION(klon,nsw), INTENT(OUT)   :: alb_dir_new,alb_dif_new
148    !albedo SB <<<     
149    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
150    REAL, DIMENSION(klon), INTENT(OUT)       :: tsurf_new    ! sea surface temperature, in K
151    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l     
152    REAL, DIMENSION(klon), INTENT(OUT)       :: lmt_bils
153    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
154
155    REAL, intent(out):: tkt(:) ! (knon)
156    ! �paisseur (m) de la couche de diffusion thermique (microlayer)
157    ! cool skin thickness
158
159    REAL, intent(out):: tks(:) ! (knon)
160    ! �paisseur (m) de la couche de diffusion de masse (microlayer)
161
162    REAL, intent(out):: taur(:) ! (knon)
163    ! momentum flux due to rain, in Pa
164
165    real, intent(out):: sss(:) ! (klon)
166    ! Bulk salinity of the surface layer of the ocean, in ppt. (Only
167    ! defined for subscripts 1:knon, but we have to declare it with
168    ! size klon because of the coupling machinery.)
169
170#ifdef ISO
171    REAL, DIMENSION(ntraciso,klon), INTENT(out) :: xtevap ! isotopes in surface evaporation flux
172    REAL, DIMENSION(klon), INTENT(out)          :: h1 ! just a diagnostic, not useful for the simulation   
173#endif
174
175    ! Local variables
176    !*************************************************************************
177    INTEGER               :: i, k
178    REAL                  :: tmp
179    REAL, PARAMETER       :: cepdu2=(0.1)**2
180    REAL, DIMENSION(klon) :: alb_eau, z0_lim
181    REAL, DIMENSION(klon) :: radsol
182    REAL, DIMENSION(klon) :: cdragq ! Cdrag pour l'evaporation
183    REAL, DIMENSION(klon) :: precip_totsnow
184    CHARACTER(len=20),PARAMETER :: modname="surf_ocean"
185    REAL rhoa(knon) ! density of moist air  (kg / m3)
186    REAL sens_prec_liq(knon)
187
188    REAL t_int(knon) ! ocean-air interface temperature, in K
189    REAL s_int(knon) ! ocean-air interface salinity, in ppt
190
191    !**************************************************************************
192
193#ifdef ISO
194#ifdef ISOVERIF
195    DO i = 1, knon
196      IF (iso_eau > 0) THEN         
197        CALL iso_verif_egalite_choix(xtspechum(iso_eau,i), &
198     &          spechum(i),'surf_ocean_mod 117', &
199     &          errmax,errmaxrel)         
200        CALL iso_verif_egalite_choix(xtsnow(iso_eau,i), &
201     &          snow(i),'surf_ocean_mod 127', &
202     &          errmax,errmaxrel)
203      ENDIF !IF (iso_eau > 0) then
204    ENDDO !DO i=1,klon
205#endif     
206#endif
207
208    !******************************************************************************
209    ! Calculate total net radiance at surface
210    !
211    !******************************************************************************
212    radsol(1:klon) = 0.0 ! initialisation a priori inutile
213    radsol(1:knon) = swnet(1:knon) + lwnet(1:knon)
214
215
216    !****************************************************************************************
217    !Total solid precip
218
219    IF (ok_bs) THEN
220       precip_totsnow(:)=precip_snow(:)+precip_bs(:)
221    ELSE
222       precip_totsnow(:)=precip_snow(:)
223    ENDIF
224
225
226    !******************************************************************************
227    ! Cdragq computed from cdrag
228    ! The difference comes only from a factor (f_z0qh_oce) on z0, so that
229    ! it can be computed inside surf_ocean
230    ! More complicated appraches may require the propagation through
231    ! pbl_surface of an independant cdragq variable.
232    !******************************************************************************
233
234    IF ( f_z0qh_oce .ne. 1.) THEN
235       ! Si on suit les formulations par exemple de Tessel, on
236       ! a z0h=0.4*nu/u*, z0q=0.62*nu/u*, d'ou f_z0qh_oce=0.62/0.4=1.55
237       cdragq(1:knon)=cdragh(1:knon)*                                      &
238            log(z1lay(1:knon)/z0h(1:knon))/log(z1lay(1:knon)/(f_z0qh_oce*z0h(1:knon)))
239    ELSE
240       cdragq(1:knon)=cdragh(1:knon)
241    ENDIF
242
243    rhoa = PS(:KNON) / (Rd * temp_air(:knon) * (1. + retv * spechum(:knon)))
244    !******************************************************************************
245    ! Switch according to type of ocean (couple, slab or forced)
246    !******************************************************************************
247    SELECT CASE(type_ocean)
248    CASE('couple')
249       CALL ocean_cpl_noice( &
250            swnet, lwnet, alb1, &
251            windsp, fder, &
252            itime, dtime, knon, knindex, &
253            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_totsnow,temp_air,spechum,&
254            AcoefH, AcoefQ, BcoefH, BcoefQ, &
255            AcoefU, AcoefV, BcoefU, BcoefV, &
256            ps, u1, v1, gustiness, tsurf_in, &
257            radsol, snow, agesno, &
258            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
259            tsurf_new, dflux_s, dflux_l, sens_prec_liq, sss, delta_sal, rhoa, &
260            delta_sst, dTer, dSer, dt_ds)
261
262    CASE('slab')
263       CALL ocean_slab_noice( &
264            itime, dtime, jour, knon, knindex, &
265            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_totsnow, temp_air, spechum,&
266            AcoefH, AcoefQ, BcoefH, BcoefQ, &
267            AcoefU, AcoefV, BcoefU, BcoefV, &
268            ps, u1, v1, gustiness, tsurf_in, &
269            radsol, snow, &
270            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
271            tsurf_new, dflux_s, dflux_l, lmt_bils)
272
273    CASE('force')
274       CALL ocean_forced_noice( &
275            itime, dtime, jour, knon, knindex, &
276            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_totsnow, &
277            temp_air, spechum, &
278            AcoefH, AcoefQ, BcoefH, BcoefQ, &
279            AcoefU, AcoefV, BcoefU, BcoefV, &
280            ps, u1, v1, gustiness, tsurf_in, &
281            radsol, snow, agesno, &
282            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
283            tsurf_new, dflux_s, dflux_l, sens_prec_liq, rhoa &
284#ifdef ISO
285            ,xtprecip_rain, xtprecip_snow, xtspechum,Roce,rlat, &
286            xtsnow,xtevap,h1 & 
287#endif           
288            )
289    END SELECT
290
291    !******************************************************************************
292    ! fcodron: compute lmt_bils  forced case (same as wfbils_oce / 1.-contfracatm)
293    !******************************************************************************
294    IF (type_ocean.NE.'slab') THEN
295       lmt_bils(1:klon)=0.
296       DO i=1,knon
297          lmt_bils(knindex(i))=(swnet(i)+lwnet(i)+fluxsens(i)+fluxlat(i)) &
298               *pctsrf(knindex(i),is_oce)/(1.-zmasq(knindex(i)))
299       END DO
300    END IF
301
302    !******************************************************************************
303    ! Calculate ocean surface albedo
304    !******************************************************************************
305    !albedo SB >>>
306    IF (iflag_albedo==0) THEN
307       !--old parametrizations of ocean surface albedo
308       !
309       IF (iflag_cycle_diurne.GE.1) THEN
310          !
311          CALL alboc_cd(rmu0,alb_eau)
312          !
313          !--ad-hoc correction for model radiative balance tuning
314          !--now outside alboc_cd routine
315          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
316          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.0),1.0)
317          !
318       ELSE
319          !
320          CALL alboc(REAL(jour),rlat,alb_eau)
321          !--ad-hoc correction for model radiative balance tuning
322          !--now outside alboc routine
323          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
324          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.04),0.60)
325          !
326       ENDIF
327       !
328       DO i =1, knon
329          DO  k=1,nsw
330             alb_dir_new(i,k) = alb_eau(knindex(i))
331          ENDDO
332       ENDDO
333       !IM 09122015 next line corresponds to the old way of doing in LMDZ5A/IPSLCM5A versions
334       !albedo for diffuse radiation is taken the same as for direct radiation
335       alb_dif_new(1:knon,:)=alb_dir_new(1:knon,:)
336       !IM 09122015 end
337       !
338    ELSE IF (iflag_albedo==1) THEN
339       !--new parametrization of ocean surface albedo by Sunghye Baek
340       !--albedo for direct and diffuse radiation are different
341       !
342       CALL ocean_albedo(knon,rmu0,knindex,windsp,SFRWL,alb_dir_new,alb_dif_new)
343       !
344       !--ad-hoc correction for model radiative balance tuning
345       alb_dir_new(1:knon,:) = fmagic*alb_dir_new(1:knon,:) + pmagic
346       alb_dif_new(1:knon,:) = fmagic*alb_dif_new(1:knon,:) + pmagic
347       alb_dir_new(1:knon,:)=MIN(MAX(alb_dir_new(1:knon,:),0.0),1.0)
348       alb_dif_new(1:knon,:)=MIN(MAX(alb_dif_new(1:knon,:),0.0),1.0)
349       !
350    ELSE IF (iflag_albedo==2) THEN
351       ! F. Codron albedo read from limit.nc
352       CALL limit_read_rug_alb(itime, dtime, jour,&
353            knon, knindex, z0_lim, alb_eau)
354       DO i =1, knon
355          DO  k=1,nsw
356             alb_dir_new(i,k) = alb_eau(i)
357          ENDDO
358       ENDDO
359       alb_dif_new=alb_dir_new
360    ENDIF
361    !albedo SB <<<
362
363    !******************************************************************************
364    ! Calculate the rugosity
365    !******************************************************************************
366    IF (iflag_z0_oce==0) THEN
367       DO i = 1, knon
368          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
369          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
370               +  0.11*14e-6 / SQRT(cdragm(i) * tmp)
371          z0m(i) = MAX(1.5e-05,z0m(i))
372       ENDDO
373       z0h(1:knon)=z0m(1:knon) ! En attendant mieux
374
375    ELSE IF (iflag_z0_oce==1) THEN
376       DO i = 1, knon
377          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
378          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
379               + 0.11*14e-6 / SQRT(cdragm(i) * tmp)
380          z0m(i) = MAX(1.5e-05,z0m(i))
381          z0h(i)=0.4*14e-6 / SQRT(cdragm(i) * tmp)
382       ENDDO
383    ELSE IF (iflag_z0_oce==-1) THEN
384       DO i = 1, knon
385          z0m(i) = z0min
386          z0h(i) = z0min
387       ENDDO
388    ELSE
389       CALL abort_physic(modname,'version non prevue',1)
390    ENDIF
391
392    if (activate_ocean_skin >= 1) then
393       if (type_ocean /= 'couple') sss(:knon) = 35.
394       call bulk_flux(tkt, tks, taur, dter, dser, t_int, s_int, ds_ns, dt_ns, &
395            u = windsp(:knon), t_ocean_1 = tsurf_new(:knon), s1 = sss(:knon), &
396            rain = precip_rain(:knon) + precip_totsnow(:knon), &
397            hf = - fluxsens(:knon), hlb = - fluxlat(:knon), &
398            rnl = - lwnet(:knon), &
399            tau = sqrt(flux_u1(:knon)**2 + flux_v1(:knon)**2), rhoa = rhoa, &
400            xlv = [(rlvtt, i = 1, knon)], rf = - sens_prec_liq, dtime = dtime, &
401            rns = swnet(:knon))
402       delta_sst = t_int - tsurf_new(:knon)
403       delta_sal = s_int - sss(:knon)
404
405       if (activate_ocean_skin == 2) then
406          tsurf_new(:knon) = t_int
407          if (type_ocean == 'couple') dt_ds = (tks / tkt) * dter
408       end if
409    end if
410   
411  END SUBROUTINE surf_ocean
412  !****************************************************************************
413  !
414END MODULE surf_ocean_mod
Note: See TracBrowser for help on using the repository browser.