source: LMDZ6/branches/Ocean_skin/libf/phylmd/ocean_forced_mod.F90 @ 3687

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

Modify sensible heat due to rain sent to the ocean

Modify the sensible heat flux due to rain which is sent to the
ocean. Replace the computation of sens_prec_liq in procedure
calcul_fluxs by a call to sens_heat_rain. Set sens_prec_sol in
procedure calcul_fluxs to 0 because, for now, sens_heat_rain is
supposed to account for both rain and snow.

For the call to sens_heat_rain in procedure calcul_fluxs, we need
an additional dummy argument rhoa of calcul_fluxs. Add dummy
argument rhoa to ocean_cpl_noice, ocean_forced_noice,
ocean_forced_ice and ocean_cpl_ice because we need to pass it down
to calcul_fluxs.

Change the dimension of sens_prec_liq and sens_prec_sol in
procedures calcul_fluxs, ocean_cpl_noice, ocean_cpl_ice,
ocean_forced_noice, ocean_forced_ice, cpl_send_ocean_fields and
cpl_send_seaice_fields to knon.

In procedures ocean_forced_noice and ocean_cpl_noice, promote
sens_prec_liq from local variable to dummy argument because we need
it in surf_ocean. Remove useless initialization of sens_prec_liq
and sens_prec_sol in ocean_cpl_noice, ocean_cpl_ice,
ocean_forced_ice and ocean_forced_noice: they are intent out in
calcul_fluxs.

Remove variable rf of module phys_output_var_mod, we use
sens_prec_liq instead. Remove local variable yrf of procedure
pbl_surface. rf and yrf appeared in pbl_surface only to be output.
Remove variable o_rf of module phys_output_ctrlout_mod. Remove
dummy argument rf of procedure surf_ocean.

Do not call sens_heat_rain in surf_ocean since we now call it
from calcul_fluxs. Move the computation of rhoa in surf_ocean
before the calls to ocean_cpl_noice and ocean_forced_noice.

Add the computation of rhoa in surf_seaice, to pass it down to
ocean_cpl_ice and ocean_forced_ice.

If activate_ocean_skin == 1 then the results are changed because the
call to sens_heat_rain in calcul_fluxs now uses the SST from the
current time step of physics. On this point, the present revision
reverses revision [3463].

  • 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 Id
File size: 13.0 KB
Line 
1!
2! $Id: ocean_forced_mod.F90 3687 2020-05-27 14:52:16Z lguez $
3!
4MODULE ocean_forced_mod
5!
6! This module is used for both the sub-surfaces ocean and sea-ice for the case of a
7! forced ocean,  "ocean=force".
8!
9  IMPLICIT NONE
10
11CONTAINS
12!
13!****************************************************************************************
14!
15  SUBROUTINE ocean_forced_noice( &
16       itime, dtime, jour, knon, knindex, &
17       p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, &
18       temp_air, spechum, &
19       AcoefH, AcoefQ, BcoefH, BcoefQ, &
20       AcoefU, AcoefV, BcoefU, BcoefV, &
21       ps, u1, v1, gustiness, tsurf_in, &
22       radsol, snow, agesno, &
23       qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
24       tsurf_new, dflux_s, dflux_l, sens_prec_liq, rhoa)
25!
26! This subroutine treats the "open ocean", all grid points that are not entierly covered
27! by ice.
28! The routine receives data from climatologie file limit.nc and does some calculations at the
29! surface.
30!
31    USE dimphy
32    USE calcul_fluxs_mod
33    USE limit_read_mod
34    USE mod_grid_phy_lmdz
35    USE indice_sol_mod
36    USE phys_output_var_mod, ONLY : sens_prec_liq_o, sens_prec_sol_o, lat_prec_liq_o, lat_prec_sol_o
37    use config_ocean_skin_m, only: activate_ocean_skin
38
39    INCLUDE "YOMCST.h"
40    INCLUDE "clesphys.h"
41
42
43! Input arguments
44!****************************************************************************************
45    INTEGER, INTENT(IN)                      :: itime, jour, knon
46    INTEGER, DIMENSION(klon), INTENT(IN)     :: knindex
47    REAL, INTENT(IN)                         :: dtime
48    REAL, DIMENSION(klon), INTENT(IN)        :: p1lay
49    REAL, DIMENSION(klon), INTENT(IN)        :: cdragh, cdragq, cdragm
50    REAL, DIMENSION(klon), INTENT(IN)        :: precip_rain, precip_snow
51    REAL, DIMENSION(klon), INTENT(IN)        :: temp_air, spechum
52    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefH, AcoefQ, BcoefH, BcoefQ
53    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefU, AcoefV, BcoefU, BcoefV
54    REAL, DIMENSION(klon), INTENT(IN)        :: ps
55    REAL, DIMENSION(klon), INTENT(IN)        :: u1, v1, gustiness
56    REAL, DIMENSION(klon), INTENT(IN)        :: tsurf_in
57    real, intent(in):: rhoa(:) ! (knon) density of moist air  (kg / m3)
58
59! In/Output arguments
60!****************************************************************************************
61    REAL, DIMENSION(klon), INTENT(INOUT)     :: radsol
62    REAL, DIMENSION(klon), INTENT(INOUT)     :: snow
63    REAL, DIMENSION(klon), INTENT(INOUT)     :: agesno !? put to 0 in ocean
64 
65! Output arguments
66!****************************************************************************************
67    REAL, DIMENSION(klon), INTENT(OUT)       :: qsurf
68    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
69    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
70    REAL, DIMENSION(klon), INTENT(OUT)       :: tsurf_new
71    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l
72    REAL, intent(out):: sens_prec_liq(:) ! (knon)
73
74! Local variables
75!****************************************************************************************
76    INTEGER                     :: i, j
77    REAL, DIMENSION(klon)       :: cal, beta, dif_grnd
78    REAL, DIMENSION(klon)       :: alb_neig, tsurf_lim, zx_sl
79    REAL, DIMENSION(klon)       :: u0, v0
80    REAL, DIMENSION(klon)       :: u1_lay, v1_lay
81    LOGICAL                     :: check=.FALSE.
82    REAL sens_prec_sol(knon)
83    REAL, DIMENSION(klon) :: lat_prec_liq, lat_prec_sol   
84
85!****************************************************************************************
86! Start calculation
87!****************************************************************************************
88    IF (check) WRITE(*,*)' Entering ocean_forced_noice'
89   
90!****************************************************************************************
91! 1)   
92! Read sea-surface temperature from file limit.nc
93!
94!****************************************************************************************
95!--sb:
96!!jyg    if (knon.eq.1) then ! single-column model
97    if (klon_glo.eq.1) then ! single-column model
98      CALL read_tsurf1d(knon,tsurf_lim) ! new
99    else ! GCM
100      CALL limit_read_sst(knon,knindex,tsurf_lim)
101    endif ! knon
102!sb--
103
104!****************************************************************************************
105! 2)
106! Flux calculation
107!
108!****************************************************************************************
109! Set some variables for calcul_fluxs
110    cal = 0.
111    beta = 1.
112    dif_grnd = 0.
113    alb_neig(:) = 0.
114    agesno(:) = 0.
115    lat_prec_liq = 0.; lat_prec_sol = 0.
116
117! Suppose zero surface speed
118    u0(:)=0.0
119    v0(:)=0.0
120    u1_lay(:) = u1(:) - u0(:)
121    v1_lay(:) = v1(:) - v0(:)
122
123! Calcul de tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l and qsurf
124    CALL calcul_fluxs(knon, is_oce, dtime, &
125         merge(tsurf_in, tsurf_lim, activate_ocean_skin == 2), p1lay, cal, &
126         beta, cdragh, cdragq, ps, &
127         precip_rain, precip_snow, snow, qsurf,  &
128         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
129         f_qsat_oce,AcoefH, AcoefQ, BcoefH, BcoefQ, &
130         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
131         sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol, rhoa)
132    if (activate_ocean_skin == 2) tsurf_new = tsurf_lim
133
134    do j = 1, knon
135      i = knindex(j)
136      sens_prec_liq_o(i,1) = sens_prec_liq(j)
137      sens_prec_sol_o(i,1) = sens_prec_sol(j)
138      lat_prec_liq_o(i,1) = lat_prec_liq(j)
139      lat_prec_sol_o(i,1) = lat_prec_sol(j)
140    enddo
141
142
143! - Flux calculation at first modele level for U and V
144    CALL calcul_flux_wind(knon, dtime, &
145         u0, v0, u1, v1, gustiness, cdragm, &
146         AcoefU, AcoefV, BcoefU, BcoefV, &
147         p1lay, temp_air, &
148         flux_u1, flux_v1) 
149
150  END SUBROUTINE ocean_forced_noice
151!
152!***************************************************************************************
153!
154  SUBROUTINE ocean_forced_ice( &
155       itime, dtime, jour, knon, knindex, &
156       tsurf_in, p1lay, cdragh, cdragm, precip_rain, precip_snow, temp_air, spechum, &
157       AcoefH, AcoefQ, BcoefH, BcoefQ, &
158       AcoefU, AcoefV, BcoefU, BcoefV, &
159       ps, u1, v1, gustiness, &
160       radsol, snow, qsol, agesno, tsoil, &
161       qsurf, alb1_new, alb2_new, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
162       tsurf_new, dflux_s, dflux_l, rhoa)
163!
164! This subroutine treats the ocean where there is ice.
165! The routine reads data from climatologie file and does flux calculations at the
166! surface.
167!
168    USE dimphy
169    USE calcul_fluxs_mod
170    USE surface_data,     ONLY : calice, calsno
171    USE limit_read_mod
172    USE fonte_neige_mod,  ONLY : fonte_neige
173    USE indice_sol_mod
174    USE phys_output_var_mod, ONLY : sens_prec_liq_o, sens_prec_sol_o, lat_prec_liq_o, lat_prec_sol_o
175
176!    INCLUDE "indicesol.h"
177    INCLUDE "dimsoil.h"
178    INCLUDE "YOMCST.h"
179    INCLUDE "clesphys.h"
180
181! Input arguments
182!****************************************************************************************
183    INTEGER, INTENT(IN)                  :: itime, jour, knon
184    INTEGER, DIMENSION(klon), INTENT(IN) :: knindex
185    REAL, INTENT(IN)                     :: dtime
186    REAL, DIMENSION(klon), INTENT(IN)    :: tsurf_in
187    REAL, DIMENSION(klon), INTENT(IN)    :: p1lay
188    REAL, DIMENSION(klon), INTENT(IN)    :: cdragh, cdragm
189    REAL, DIMENSION(klon), INTENT(IN)    :: precip_rain, precip_snow
190    REAL, DIMENSION(klon), INTENT(IN)    :: temp_air, spechum
191    REAL, DIMENSION(klon), INTENT(IN)    :: AcoefH, AcoefQ, BcoefH, BcoefQ
192    REAL, DIMENSION(klon), INTENT(IN)    :: AcoefU, AcoefV, BcoefU, BcoefV
193    REAL, DIMENSION(klon), INTENT(IN)    :: ps
194    REAL, DIMENSION(klon), INTENT(IN)    :: u1, v1, gustiness
195    real, intent(in):: rhoa(:) ! (knon) density of moist air  (kg / m3)
196
197! In/Output arguments
198!****************************************************************************************
199    REAL, DIMENSION(klon), INTENT(INOUT)          :: radsol
200    REAL, DIMENSION(klon), INTENT(INOUT)          :: snow, qsol
201    REAL, DIMENSION(klon), INTENT(INOUT)          :: agesno
202    REAL, DIMENSION(klon, nsoilmx), INTENT(INOUT) :: tsoil
203
204! Output arguments
205!****************************************************************************************
206    REAL, DIMENSION(klon), INTENT(OUT)            :: qsurf
207    REAL, DIMENSION(klon), INTENT(OUT)            :: alb1_new  ! new albedo in visible SW interval
208    REAL, DIMENSION(klon), INTENT(OUT)            :: alb2_new  ! new albedo in near IR interval
209    REAL, DIMENSION(klon), INTENT(OUT)            :: evap, fluxsens, fluxlat
210    REAL, DIMENSION(klon), INTENT(OUT)            :: flux_u1, flux_v1
211    REAL, DIMENSION(klon), INTENT(OUT)            :: tsurf_new
212    REAL, DIMENSION(klon), INTENT(OUT)            :: dflux_s, dflux_l     
213
214! Local variables
215!****************************************************************************************
216    LOGICAL                     :: check=.FALSE.
217    INTEGER                     :: i, j
218    REAL                        :: zfra
219    REAL, PARAMETER             :: t_grnd=271.35
220    REAL, DIMENSION(klon)       :: cal, beta, dif_grnd, capsol
221    REAL, DIMENSION(klon)       :: alb_neig, tsurf_tmp
222    REAL, DIMENSION(klon)       :: soilcap, soilflux
223    REAL, DIMENSION(klon)       :: u0, v0
224    REAL, DIMENSION(klon)       :: u1_lay, v1_lay
225    REAL sens_prec_liq(knon), sens_prec_sol (knon)   
226    REAL, DIMENSION(klon)       :: lat_prec_liq, lat_prec_sol   
227
228
229!****************************************************************************************
230! Start calculation
231!****************************************************************************************
232    IF (check) WRITE(*,*)'Entering surface_seaice, knon=',knon
233
234!****************************************************************************************
235! 1)
236! Flux calculation : tsurf_new, evap, fluxlat, fluxsens, flux_u1, flux_v1
237!                    dflux_s, dflux_l and qsurf
238!****************************************************************************************
239
240    tsurf_tmp(:) = tsurf_in(:)
241
242! calculate the parameters cal, beta, capsol and dif_grnd
243    CALL calbeta(dtime, is_sic, knon, snow, qsol, beta, capsol, dif_grnd)
244
245   
246    IF (soil_model) THEN
247! update tsoil and calculate soilcap and soilflux
248       CALL soil(dtime, is_sic, knon, snow, tsurf_tmp, tsoil,soilcap, soilflux)
249       cal(1:knon) = RCPD / soilcap(1:knon)
250       radsol(1:knon) = radsol(1:knon)  + soilflux(1:knon)
251       dif_grnd = 1.0 / tau_gl
252    ELSE
253       dif_grnd = 1.0 / tau_gl
254       cal = RCPD * calice
255       WHERE (snow > 0.0) cal = RCPD * calsno
256    ENDIF
257
258    beta = 1.0
259    lat_prec_liq = 0.; lat_prec_sol = 0.
260
261! Suppose zero surface speed
262    u0(:)=0.0
263    v0(:)=0.0
264    u1_lay(:) = u1(:) - u0(:)
265    v1_lay(:) = v1(:) - v0(:)
266    CALL calcul_fluxs(knon, is_sic, dtime, &
267         tsurf_tmp, p1lay, cal, beta, cdragh, cdragh, ps, &
268         precip_rain, precip_snow, snow, qsurf,  &
269         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
270         f_qsat_oce,AcoefH, AcoefQ, BcoefH, BcoefQ, &
271         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
272         sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol, rhoa)
273    do j = 1, knon
274      i = knindex(j)
275      sens_prec_liq_o(i,2) = sens_prec_liq(j)
276      sens_prec_sol_o(i,2) = sens_prec_sol(j)
277      lat_prec_liq_o(i,2) = lat_prec_liq(j)
278      lat_prec_sol_o(i,2) = lat_prec_sol(j)
279    enddo
280
281! - Flux calculation at first modele level for U and V
282    CALL calcul_flux_wind(knon, dtime, &
283         u0, v0, u1, v1, gustiness, cdragm, &
284         AcoefU, AcoefV, BcoefU, BcoefV, &
285         p1lay, temp_air, &
286         flux_u1, flux_v1) 
287
288!****************************************************************************************
289! 2)
290! Calculations due to snow and runoff
291!
292!****************************************************************************************
293    CALL fonte_neige( knon, is_sic, knindex, dtime, &
294         tsurf_tmp, precip_rain, precip_snow, &
295         snow, qsol, tsurf_new, evap)
296   
297! Calculation of albedo at snow (alb_neig) and update the age of snow (agesno)
298!
299    CALL albsno(klon, knon, dtime, agesno(:), alb_neig(:), precip_snow(:)) 
300
301    WHERE (snow(1:knon) .LT. 0.0001) agesno(1:knon) = 0.
302
303    alb1_new(:) = 0.0
304    DO i=1, knon
305       zfra = MAX(0.0,MIN(1.0,snow(i)/(snow(i)+10.0)))
306       alb1_new(i) = alb_neig(i) * zfra +  0.6 * (1.0-zfra)
307    ENDDO
308
309    alb2_new(:) = alb1_new(:)
310
311  END SUBROUTINE ocean_forced_ice
312
313!************************************************************************
314! 1D case
315!************************************************************************
316  SUBROUTINE read_tsurf1d(knon,sst_out)
317
318! This subroutine specifies the surface temperature to be used in 1D simulations
319
320      USE dimphy, ONLY : klon
321
322      INTEGER, INTENT(IN)                  :: knon     ! nomber of points on compressed grid
323      REAL, DIMENSION(klon), INTENT(OUT)   :: sst_out  ! tsurf used to force the single-column model
324
325       INTEGER :: i
326! COMMON defined in lmdz1d.F:
327       real ts_cur
328       common /sst_forcing/ts_cur
329
330       DO i = 1, knon
331        sst_out(i) = ts_cur
332       ENDDO
333
334      END SUBROUTINE read_tsurf1d
335
336!
337!************************************************************************
338!
339END MODULE ocean_forced_mod
340
341
342
343
344
345
Note: See TracBrowser for help on using the repository browser.