GCC Code Coverage Report
Directory: ./ Exec Total Coverage
File: phylmd/surf_ocean_mod.F90 Lines: 41 80 51.2 %
Date: 2023-06-30 12:51:15 Branches: 52 152 34.2 %

Line Branch Exec Source
1
!
2
! $Id: surf_ocean_mod.F90 4526 2023-05-08 12:35:08Z evignon $
3
!
4
MODULE surf_ocean_mod
5
6
  IMPLICIT NONE
7
8
CONTAINS
9
  !
10
  !******************************************************************************
11
  !
12
219456
  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
288
       ps, u1, v1, gustiness, rugoro, pctsrf, &
19
       snow, qsurf, agesno, &
20
288
       z0m, z0h, SFRWL, alb_dir_new, alb_dif_new, evap, fluxsens, fluxlat, &
21
       tsurf_new, dflux_s, dflux_l, lmt_bils, &
22


576
       flux_u1, flux_v1, delta_sst, delta_sal, ds_ns, dt_ns, dter, dser, &
23



576
       dt_ds, tkt, tks, 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, DIMENSION(klon), INTENT(IN)        :: tsurf_in     ! 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, precip_bs
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. Defined
84
    ! only if activate_ocean_skin >= 1.
85
86
    REAL, intent(inout):: ds_ns(:) ! (knon)
87
    ! "delta salinity near surface". Salinity variation in the
88
    ! near-surface turbulent layer. That is subskin salinity minus
89
    ! foundation salinity. In ppt.
90
91
    REAL, intent(inout):: dt_ns(:) ! (knon)
92
    ! "delta temperature near surface". Temperature variation in the
93
    ! near-surface turbulent layer. That is subskin temperature
94
    ! minus foundation temperature. (Can be negative.) In K.
95
96
    REAL, intent(inout):: dter(:) ! (knon)
97
    ! Temperature variation in the diffusive microlayer, that is
98
    ! ocean-air interface temperature minus subskin temperature. In
99
    ! K.
100
101
    REAL, intent(inout):: dser(:) ! (knon)
102
    ! Salinity variation in the diffusive microlayer, that is
103
    ! ocean-air interface salinity minus subskin salinity. In ppt.
104
105
    real, intent(inout):: dt_ds(:) ! (knon)
106
    ! (tks / tkt) * dTer, in K
107
108
    ! Output variables
109
    !**************************************************************************
110
    REAL, DIMENSION(klon), INTENT(OUT)       :: z0m
111
    !albedo SB >>>
112
    !    REAL, DIMENSION(klon), INTENT(OUT)  :: alb1_new  ! new albedo in visible SW interval
113
    !    REAL, DIMENSION(klon), INTENT(OUT)  :: alb2_new  ! new albedo in near IR interval
114
    REAL, DIMENSION(6), INTENT(IN)           :: SFRWL
115
    REAL, DIMENSION(klon,nsw), INTENT(OUT)   :: alb_dir_new,alb_dif_new
116
    !albedo SB <<<
117
    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
118
    REAL, DIMENSION(klon), INTENT(OUT)       :: tsurf_new    ! sea surface temperature, in K
119
    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l
120
    REAL, DIMENSION(klon), INTENT(OUT)       :: lmt_bils
121
    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
122
123
    REAL, intent(out):: tkt(:) ! (knon)
124
    ! �paisseur (m) de la couche de diffusion thermique (microlayer)
125
    ! cool skin thickness
126
127
    REAL, intent(out):: tks(:) ! (knon)
128
    ! �paisseur (m) de la couche de diffusion de masse (microlayer)
129
130
    REAL, intent(out):: taur(:) ! (knon)
131
    ! momentum flux due to rain, in Pa
132
133
    real, intent(out):: sss(:) ! (klon)
134
    ! Bulk salinity of the surface layer of the ocean, in ppt. (Only
135
    ! defined for subscripts 1:knon, but we have to declare it with
136
    ! size klon because of the coupling machinery.)
137
138
    ! Local variables
139
    !*************************************************************************
140
    INTEGER               :: i, k
141
    REAL                  :: tmp
142
    REAL, PARAMETER       :: cepdu2=(0.1)**2
143
576
    REAL, DIMENSION(klon) :: alb_eau, z0_lim
144
576
    REAL, DIMENSION(klon) :: radsol
145
576
    REAL, DIMENSION(klon) :: cdragq ! Cdrag pour l'evaporation
146
576
    REAL, DIMENSION(klon) :: precip_totsnow
147
    CHARACTER(len=20),PARAMETER :: modname="surf_ocean"
148
576
    real rhoa(knon) ! density of moist air  (kg / m3)
149
576
    REAL sens_prec_liq(knon)
150
151
576
    REAL t_int(knon) ! ocean-air interface temperature, in K
152
288
    real s_int(knon) ! ocean-air interface salinity, in ppt
153
154
    !**************************************************************************
155
156
157
    !******************************************************************************
158
    ! Calculate total net radiance at surface
159
    !
160
    !******************************************************************************
161
286560
    radsol(1:klon) = 0.0 ! initialisation a priori inutile
162
218016
    radsol(1:knon) = swnet(1:knon) + lwnet(1:knon)
163
164
165
    !****************************************************************************************
166
    !Total solid precip
167
168
288
    IF (ok_bs) THEN
169
       precip_totsnow(:)=precip_snow(:)+precip_bs(:)
170
    ELSE
171
286560
       precip_totsnow(:)=precip_snow(:)
172
    ENDIF
173
174
175
    !******************************************************************************
176
    ! Cdragq computed from cdrag
177
    ! The difference comes only from a factor (f_z0qh_oce) on z0, so that
178
    ! it can be computed inside surf_ocean
179
    ! More complicated appraches may require the propagation through
180
    ! pbl_surface of an independant cdragq variable.
181
    !******************************************************************************
182
183
288
    IF ( f_z0qh_oce .ne. 1.) THEN
184
       ! Si on suit les formulations par exemple de Tessel, on
185
       ! a z0h=0.4*nu/u*, z0q=0.62*nu/u*, d'ou f_z0qh_oce=0.62/0.4=1.55
186
       cdragq(1:knon)=cdragh(1:knon)*                                      &
187
218016
            log(z1lay(1:knon)/z0h(1:knon))/log(z1lay(1:knon)/(f_z0qh_oce*z0h(1:knon)))
188
    ELSE
189
       cdragq(1:knon)=cdragh(1:knon)
190
    ENDIF
191
192
218016
    rhoa = PS(:KNON) / (Rd * temp_air(:knon) * (1. + retv * spechum(:knon)))
193
    !******************************************************************************
194
    ! Switch according to type of ocean (couple, slab or forced)
195
    !******************************************************************************
196
    SELECT CASE(type_ocean)
197
    CASE('couple')
198
       CALL ocean_cpl_noice( &
199
            swnet, lwnet, alb1, &
200
            windsp, fder, &
201
            itime, dtime, knon, knindex, &
202
            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_totsnow,temp_air,spechum,&
203
            AcoefH, AcoefQ, BcoefH, BcoefQ, &
204
            AcoefU, AcoefV, BcoefU, BcoefV, &
205
            ps, u1, v1, gustiness, tsurf_in, &
206
            radsol, snow, agesno, &
207
            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
208
            tsurf_new, dflux_s, dflux_l, sens_prec_liq, sss, delta_sal, rhoa, &
209
            delta_sst, dTer, dSer, dt_ds)
210
211
    CASE('slab')
212
       CALL ocean_slab_noice( &
213
            itime, dtime, jour, knon, knindex, &
214
            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_totsnow, temp_air, spechum,&
215
            AcoefH, AcoefQ, BcoefH, BcoefQ, &
216
            AcoefU, AcoefV, BcoefU, BcoefV, &
217
            ps, u1, v1, gustiness, tsurf_in, &
218
            radsol, snow, &
219
            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
220
            tsurf_new, dflux_s, dflux_l, lmt_bils)
221
222
    CASE('force')
223
       CALL ocean_forced_noice( &
224
            itime, dtime, jour, knon, knindex, &
225
            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_totsnow, &
226
            temp_air, spechum, &
227
            AcoefH, AcoefQ, BcoefH, BcoefQ, &
228
            AcoefU, AcoefV, BcoefU, BcoefV, &
229
            ps, u1, v1, gustiness, tsurf_in, &
230
            radsol, snow, agesno, &
231
            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
232

288
            tsurf_new, dflux_s, dflux_l, sens_prec_liq, rhoa)
233
    END SELECT
234
235
    !******************************************************************************
236
    ! fcodron: compute lmt_bils  forced case (same as wfbils_oce / 1.-contfracatm)
237
    !******************************************************************************
238
288
    IF (type_ocean.NE.'slab') THEN
239
286560
       lmt_bils(1:klon)=0.
240
218016
       DO i=1,knon
241
          lmt_bils(knindex(i))=(swnet(i)+lwnet(i)+fluxsens(i)+fluxlat(i)) &
242
218016
               *pctsrf(knindex(i),is_oce)/(1.-zmasq(knindex(i)))
243
       END DO
244
    END IF
245
246
    !******************************************************************************
247
    ! Calculate ocean surface albedo
248
    !******************************************************************************
249
    !albedo SB >>>
250
288
    IF (iflag_albedo==0) THEN
251
       !--old parametrizations of ocean surface albedo
252
       !
253
       IF (iflag_cycle_diurne.GE.1) THEN
254
          !
255
          CALL alboc_cd(rmu0,alb_eau)
256
          !
257
          !--ad-hoc correction for model radiative balance tuning
258
          !--now outside alboc_cd routine
259
          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
260
          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.0),1.0)
261
          !
262
       ELSE
263
          !
264
          CALL alboc(REAL(jour),rlat,alb_eau)
265
          !--ad-hoc correction for model radiative balance tuning
266
          !--now outside alboc routine
267
          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
268
          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.04),0.60)
269
          !
270
       ENDIF
271
       !
272
       DO i =1, knon
273
          DO  k=1,nsw
274
             alb_dir_new(i,k) = alb_eau(knindex(i))
275
          ENDDO
276
       ENDDO
277
       !IM 09122015 next line corresponds to the old way of doing in LMDZ5A/IPSLCM5A versions
278
       !albedo for diffuse radiation is taken the same as for direct radiation
279
       alb_dif_new(1:knon,:)=alb_dir_new(1:knon,:)
280
       !IM 09122015 end
281
       !
282
288
    ELSE IF (iflag_albedo==1) THEN
283
       !--new parametrization of ocean surface albedo by Sunghye Baek
284
       !--albedo for direct and diffuse radiation are different
285
       !
286
288
       CALL ocean_albedo(knon,rmu0,knindex,windsp,SFRWL,alb_dir_new,alb_dif_new)
287
       !
288
       !--ad-hoc correction for model radiative balance tuning
289

1308384
       alb_dir_new(1:knon,:) = fmagic*alb_dir_new(1:knon,:) + pmagic
290

1308384
       alb_dif_new(1:knon,:) = fmagic*alb_dif_new(1:knon,:) + pmagic
291

1308384
       alb_dir_new(1:knon,:)=MIN(MAX(alb_dir_new(1:knon,:),0.0),1.0)
292

1308384
       alb_dif_new(1:knon,:)=MIN(MAX(alb_dif_new(1:knon,:),0.0),1.0)
293
       !
294
    ELSE IF (iflag_albedo==2) THEN
295
       ! F. Codron albedo read from limit.nc
296
       CALL limit_read_rug_alb(itime, dtime, jour,&
297
            knon, knindex, z0_lim, alb_eau)
298
       DO i =1, knon
299
          DO  k=1,nsw
300
             alb_dir_new(i,k) = alb_eau(i)
301
          ENDDO
302
       ENDDO
303
       alb_dif_new=alb_dir_new
304
    ENDIF
305
    !albedo SB <<<
306
307
    !******************************************************************************
308
    ! Calculate the rugosity
309
    !******************************************************************************
310
288
    IF (iflag_z0_oce==0) THEN
311
       DO i = 1, knon
312
          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
313
          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
314
               +  0.11*14e-6 / SQRT(cdragm(i) * tmp)
315
          z0m(i) = MAX(1.5e-05,z0m(i))
316
       ENDDO
317
       z0h(1:knon)=z0m(1:knon) ! En attendant mieux
318
319
288
    ELSE IF (iflag_z0_oce==1) THEN
320
218016
       DO i = 1, knon
321
217728
          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
322
          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
323
217728
               + 0.11*14e-6 / SQRT(cdragm(i) * tmp)
324
217728
          z0m(i) = MAX(1.5e-05,z0m(i))
325
218016
          z0h(i)=0.4*14e-6 / SQRT(cdragm(i) * tmp)
326
       ENDDO
327
    ELSE IF (iflag_z0_oce==-1) THEN
328
       DO i = 1, knon
329
          z0m(i) = z0min
330
          z0h(i) = z0min
331
       ENDDO
332
    ELSE
333
       CALL abort_physic(modname,'version non prevue',1)
334
    ENDIF
335
336
288
    if (activate_ocean_skin >= 1) then
337
       if (type_ocean /= 'couple') sss(:knon) = 35.
338
       call bulk_flux(tkt, tks, taur, dter, dser, t_int, s_int, ds_ns, dt_ns, &
339
            u = windsp(:knon), t_ocean_1 = tsurf_new(:knon), s1 = sss(:knon), &
340
            rain = precip_rain(:knon) + precip_totsnow(:knon), &
341
            hf = - fluxsens(:knon), hlb = - fluxlat(:knon), &
342
            rnl = - lwnet(:knon), &
343
            tau = sqrt(flux_u1(:knon)**2 + flux_v1(:knon)**2), rhoa = rhoa, &
344
            xlv = [(rlvtt, i = 1, knon)], rf = - sens_prec_liq, dtime = dtime, &
345
            rns = swnet(:knon))
346
       delta_sst = t_int - tsurf_new(:knon)
347
       delta_sal = s_int - sss(:knon)
348
349
       if (activate_ocean_skin == 2) then
350
          tsurf_new(:knon) = t_int
351
          if (type_ocean == 'couple') dt_ds = (tks / tkt) * dter
352
       end if
353
    end if
354
355
288
  END SUBROUTINE surf_ocean
356
  !****************************************************************************
357
  !
358
END MODULE surf_ocean_mod