source: LMDZ6/branches/Ocean_skin/libf/phylmd/ocean_cpl_mod.F90

Last change on this file was 4020, checked in by lguez, 2 years ago

Send 3 more fields to the ocean

Send 3 more fields to the ocean to compute CO2 flux at
ocean-atmosphere interface. The three fields are dter and dser, which
already existed, and a newly created field: dt_ds. So dter and dser
have to become state variables. The variable dt_ds of module
phys_state_var_mod is only allocated and defined if
activate_ocean_skin == 2 and type_ocean == "couple".

  • 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: 15.6 KB
Line 
1!
2! $Id: ocean_cpl_mod.F90 4020 2021-11-26 07:27:26Z asima $
3!
4MODULE ocean_cpl_mod
5!
6! This module is used both for the sub-surface ocean and sea-ice for the case of a
7! coupled model configuration, ocean=couple.
8!
9
10  IMPLICIT NONE
11  PRIVATE
12
13  PUBLIC :: ocean_cpl_init, ocean_cpl_noice, ocean_cpl_ice
14
15
16!****************************************************************************************
17!
18CONTAINS
19!
20!****************************************************************************************
21!
22  SUBROUTINE ocean_cpl_init(dtime, rlon, rlat)
23!
24! Allocate fields for this module and initailize the module mod_cpl
25!
26    USE dimphy,           ONLY : klon
27    USE cpl_mod
28
29! Input arguments
30!*************************************************************************************
31    REAL, INTENT(IN)                  :: dtime
32    REAL, DIMENSION(klon), INTENT(IN) :: rlon, rlat
33
34! Local variables
35!*************************************************************************************
36    INTEGER              :: error
37    CHARACTER (len = 80) :: abort_message
38    CHARACTER (len = 20) :: modname = 'ocean_cpl_init'
39
40! Initialize module cpl_init
41    CALL cpl_init(dtime, rlon, rlat)
42   
43  END SUBROUTINE ocean_cpl_init
44!
45!****************************************************************************************
46!
47  SUBROUTINE ocean_cpl_noice( &
48       swnet, lwnet, alb1, &
49       windsp, fder_old, &
50       itime, dtime, knon, knindex, &
51       p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, temp_air, spechum, &
52       AcoefH, AcoefQ, BcoefH, BcoefQ, &
53       AcoefU, AcoefV, BcoefU, BcoefV, &
54       ps, u1, v1, gustiness, tsurf_in, &
55       radsol, snow, agesno, &
56       qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
57       tsurf_new, dflux_s, dflux_l, sens_prec_liq, sss, delta_sal, rhoa, &
58       delta_sst, dTer, dSer, dt_ds)
59
60!
61! This subroutine treats the "open ocean", all grid points that are not entierly covered
62! by ice. The subroutine first receives fields from coupler, then some calculations at
63! surface is done and finally it sends some fields to the coupler.
64!
65    USE dimphy,           ONLY : klon
66    USE calcul_fluxs_mod
67    USE indice_sol_mod
68    USE phys_output_var_mod, ONLY : sens_prec_liq_o, sens_prec_sol_o, lat_prec_liq_o, lat_prec_sol_o
69    USE cpl_mod, ONLY : gath2cpl, cpl_receive_ocean_fields, &
70         cpl_send_ocean_fields
71    use config_ocean_skin_m, only: activate_ocean_skin
72
73    INCLUDE "YOMCST.h"
74    INCLUDE "clesphys.h"
75!   
76! Input arguments 
77!****************************************************************************************
78    INTEGER, INTENT(IN)                      :: itime, knon
79    INTEGER, DIMENSION(klon), INTENT(IN)     :: knindex
80    REAL, INTENT(IN)                         :: dtime
81    REAL, DIMENSION(klon), INTENT(IN)        :: swnet
82    REAL, DIMENSION(klon), INTENT(IN)        :: lwnet
83    REAL, DIMENSION(klon), INTENT(IN)        :: alb1   ! albedo in visible SW interval
84    REAL, DIMENSION(klon), INTENT(IN)        :: windsp
85    REAL, DIMENSION(klon), INTENT(IN)        :: fder_old
86    REAL, DIMENSION(klon), INTENT(IN)        :: p1lay
87    REAL, DIMENSION(klon), INTENT(IN)        :: cdragh, cdragq, cdragm
88    REAL, DIMENSION(klon), INTENT(IN)        :: precip_rain, precip_snow
89    REAL, DIMENSION(klon), INTENT(IN)        :: temp_air, spechum
90    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefH, AcoefQ, BcoefH, BcoefQ
91    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefU, AcoefV, BcoefU, BcoefV
92    REAL, DIMENSION(klon), INTENT(IN)        :: ps
93    REAL, DIMENSION(klon), INTENT(IN)        :: u1, v1, gustiness
94    REAL, INTENT(IN) :: tsurf_in(:) ! (klon)
95
96    real, intent(in):: delta_sal(:) ! (knon)
97    ! ocean-air interface salinity minus bulk salinity, in ppt
98
99    real, intent(in):: rhoa(:) ! (knon) density of moist air  (kg / m3)
100
101    REAL, intent(in):: delta_sst(:) ! (knon)
102    ! Ocean-air interface temperature minus bulk SST, in K. Defined
103    ! only if activate_ocean_skin >= 1.
104
105    REAL, intent(in):: dter(:) ! (knon)
106    ! Temperature variation in the diffusive microlayer, that is
107    ! ocean-air interface temperature minus subskin temperature. In
108    ! K.
109
110    REAL, intent(in):: dser(:) ! (knon)
111    ! Salinity variation in the diffusive microlayer, that is
112    ! ocean-air interface salinity minus subskin salinity. In ppt.
113
114    real, intent(in):: dt_ds(:) ! (knon)
115    ! (tks / tkt) * dTer, in K
116
117! In/Output arguments
118!****************************************************************************************
119    REAL, DIMENSION(klon), INTENT(INOUT)     :: radsol
120    REAL, DIMENSION(klon), INTENT(INOUT)     :: snow
121    REAL, DIMENSION(klon), INTENT(INOUT)     :: agesno
122 
123! Output arguments
124!****************************************************************************************
125    REAL, DIMENSION(klon), INTENT(OUT)       :: qsurf
126    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
127    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
128    REAL, DIMENSION(klon), INTENT(OUT)       :: tsurf_new
129    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l     
130    REAL, intent(out):: sens_prec_liq(:) ! (knon)
131
132    REAL, INTENT(OUT):: sss(:) ! (klon)
133    ! bulk salinity of the surface layer of the ocean, in ppt
134 
135
136! Local variables
137!****************************************************************************************
138    INTEGER               :: i, j
139    INTEGER, DIMENSION(1) :: iloc
140    REAL, DIMENSION(klon) :: cal, beta, dif_grnd
141    REAL, DIMENSION(klon) :: fder_new
142    REAL, DIMENSION(klon) :: tsurf_cpl
143    REAL, DIMENSION(klon) :: u0_cpl, v0_cpl
144    REAL, DIMENSION(klon) :: u1_lay, v1_lay
145    LOGICAL               :: check=.FALSE.
146    REAL sens_prec_sol(knon) 
147    REAL, DIMENSION(klon) :: lat_prec_liq, lat_prec_sol   
148
149! End definitions
150!****************************************************************************************
151
152    IF (check) WRITE(*,*)' Entering ocean_cpl_noice'
153
154!****************************************************************************************
155! Receive sea-surface temperature(tsurf_cpl) from coupler
156!
157!****************************************************************************************
158    CALL cpl_receive_ocean_fields(knon, knindex, tsurf_cpl, u0_cpl, v0_cpl, &
159         sss)
160
161!****************************************************************************************
162! Calculate fluxes at surface
163!
164!****************************************************************************************
165    cal = 0.
166    beta = 1.
167    dif_grnd = 0.
168    agesno(:) = 0.
169    lat_prec_liq = 0.; lat_prec_sol = 0.
170   
171
172    DO i = 1, knon
173       u1_lay(i) = u1(i) - u0_cpl(i)
174       v1_lay(i) = v1(i) - v0_cpl(i)
175    END DO
176
177    CALL calcul_fluxs(knon, is_oce, dtime, &
178         merge(tsurf_in, tsurf_cpl, activate_ocean_skin == 2), p1lay, cal, &
179         beta, cdragh, cdragq, ps, &
180         precip_rain, precip_snow, snow, qsurf,  &
181         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
182         f_qsat_oce,AcoefH, AcoefQ, BcoefH, BcoefQ, &
183         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
184         sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol, rhoa)
185
186    if (activate_ocean_skin == 2) then
187       ! tsurf_new was set to tsurf_in in calcul_flux, correct it to
188       ! the new bulk SST tsurf_cpl:
189       tsurf_new = tsurf_cpl
190    end if
191
192    ! assertion: tsurf_new == tsurf_cpl
193   
194    do j = 1, knon
195      i = knindex(j)
196      sens_prec_liq_o(i,1) = sens_prec_liq(j)
197      sens_prec_sol_o(i,1) = sens_prec_sol(j)
198      lat_prec_liq_o(i,1) = lat_prec_liq(j)
199      lat_prec_sol_o(i,1) = lat_prec_sol(j)
200    enddo
201
202
203   
204! - Flux calculation at first modele level for U and V
205    CALL calcul_flux_wind(knon, dtime, &
206         u0_cpl, v0_cpl, u1, v1, gustiness, cdragm, &
207         AcoefU, AcoefV, BcoefU, BcoefV, &
208         p1lay, temp_air, &
209         flux_u1, flux_v1) 
210
211!****************************************************************************************
212! Calculate fder : flux derivative (sensible and latente)
213!
214!****************************************************************************************
215    fder_new(:) = fder_old(:) + dflux_s(:) + dflux_l(:)
216   
217    iloc = MAXLOC(fder_new(1:klon))
218    IF (check .AND. fder_new(iloc(1))> 0.) THEN
219       WRITE(*,*)'**** Debug fder****'
220       WRITE(*,*)'max fder(',iloc(1),') = ',fder_new(iloc(1))
221       WRITE(*,*)'fder_old, dflux_s, dflux_l',fder_old(iloc(1)), &
222            dflux_s(iloc(1)), dflux_l(iloc(1))
223    ENDIF
224
225!****************************************************************************************
226! Send and cumulate fields to the coupler
227!
228!****************************************************************************************
229
230    CALL cpl_send_ocean_fields(itime, knon, knindex, swnet, lwnet, fluxlat, &
231         fluxsens, precip_rain, precip_snow, evap, tsurf_new, fder_new, alb1, &
232         flux_u1, flux_v1, windsp, sens_prec_liq, sens_prec_sol, lat_prec_liq, &
233         lat_prec_sol, delta_sst, delta_sal, dTer, dSer, dt_ds)
234
235  END SUBROUTINE ocean_cpl_noice
236!
237!****************************************************************************************
238!
239  SUBROUTINE ocean_cpl_ice( &
240       rlon, rlat, swnet, lwnet, alb1, &
241       fder_old, &
242       itime, dtime, knon, knindex, &
243       lafin, &
244       p1lay, cdragh, cdragm, precip_rain, precip_snow, temp_air, spechum, &
245       AcoefH, AcoefQ, BcoefH, BcoefQ, &
246       AcoefU, AcoefV, BcoefU, BcoefV, &
247       ps, u1, v1, gustiness, pctsrf, &
248       radsol, snow, qsurf, &
249       alb1_new, alb2_new, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
250       tsurf_new, dflux_s, dflux_l, rhoa)
251!
252! This subroutine treats the ocean where there is ice. The subroutine first receives
253! fields from coupler, then some calculations at surface is done and finally sends
254! some fields to the coupler.
255!   
256    USE dimphy,           ONLY : klon
257    USE cpl_mod
258    USE calcul_fluxs_mod
259    USE indice_sol_mod
260    USE phys_output_var_mod, ONLY : sens_prec_liq_o, sens_prec_sol_o, lat_prec_liq_o, lat_prec_sol_o
261
262    INCLUDE "YOMCST.h"
263    INCLUDE "clesphys.h"
264
265! Input arguments
266!****************************************************************************************
267    INTEGER, INTENT(IN)                      :: itime, knon
268    INTEGER, DIMENSION(klon), INTENT(IN)     :: knindex
269    LOGICAL, INTENT(IN)                      :: lafin
270    REAL, INTENT(IN)                         :: dtime
271    REAL, DIMENSION(klon), INTENT(IN)        :: rlon, rlat
272    REAL, DIMENSION(klon), INTENT(IN)        :: swnet
273    REAL, DIMENSION(klon), INTENT(IN)        :: lwnet
274    REAL, DIMENSION(klon), INTENT(IN)        :: alb1   ! albedo in visible SW interval
275    REAL, DIMENSION(klon), INTENT(IN)        :: fder_old
276    REAL, DIMENSION(klon), INTENT(IN)        :: p1lay
277    REAL, DIMENSION(klon), INTENT(IN)        :: cdragh, cdragm
278    REAL, DIMENSION(klon), INTENT(IN)        :: precip_rain, precip_snow
279    REAL, DIMENSION(klon), INTENT(IN)        :: temp_air, spechum
280    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefH, AcoefQ, BcoefH, BcoefQ
281    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefU, AcoefV, BcoefU, BcoefV
282    REAL, DIMENSION(klon), INTENT(IN)        :: ps
283    REAL, DIMENSION(klon), INTENT(IN)        :: u1, v1, gustiness
284    REAL, DIMENSION(klon,nbsrf), INTENT(IN)  :: pctsrf
285    real, intent(in):: rhoa(:) ! (knon) density of moist air  (kg / m3)
286
287! In/output arguments
288!****************************************************************************************
289    REAL, DIMENSION(klon), INTENT(INOUT)     :: radsol
290    REAL, DIMENSION(klon), INTENT(INOUT)     :: snow
291
292! Output arguments
293!****************************************************************************************
294    REAL, DIMENSION(klon), INTENT(OUT)       :: qsurf
295    REAL, DIMENSION(klon), INTENT(OUT)       :: alb1_new, alb2_new
296    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
297    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
298    REAL, DIMENSION(klon), INTENT(OUT)       :: tsurf_new
299    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l     
300 
301
302! Local variables
303!****************************************************************************************
304    INTEGER                 :: i, j
305    INTEGER, DIMENSION(1)   :: iloc
306    LOGICAL                 :: check=.FALSE.
307    REAL, PARAMETER         :: t_grnd=271.35
308    REAL, DIMENSION(klon)   :: cal, beta, dif_grnd
309    REAL, DIMENSION(klon)   :: tsurf_cpl, fder_new
310    REAL, DIMENSION(klon)   :: alb_cpl
311    REAL, DIMENSION(klon)   :: u0, v0
312    REAL, DIMENSION(klon)   :: u1_lay, v1_lay
313    REAL sens_prec_liq(knon), sens_prec_sol(knon)   
314    REAL, DIMENSION(klon) :: lat_prec_liq, lat_prec_sol   
315
316! End definitions
317!****************************************************************************************
318   
319    IF (check) WRITE(*,*)'Entering surface_seaice, knon=',knon
320
321    lat_prec_liq = 0.; lat_prec_sol = 0.
322
323!****************************************************************************************
324! Receive ocean temperature(tsurf_cpl) and albedo(alb_new) from coupler
325!
326!****************************************************************************************
327
328    CALL cpl_receive_seaice_fields(knon, knindex, &
329         tsurf_cpl, alb_cpl, u0, v0)
330
331    alb1_new(1:knon) = alb_cpl(1:knon)
332    alb2_new(1:knon) = alb_cpl(1:knon)   
333
334   
335!****************************************************************************************
336! Calculate fluxes at surface
337!
338!****************************************************************************************
339    cal = 0.
340    dif_grnd = 0.
341    beta = 1.0
342   
343    DO i = 1, knon
344       u1_lay(i) = u1(i) - u0(i)
345       v1_lay(i) = v1(i) - v0(i)
346    END DO
347
348    CALL calcul_fluxs(knon, is_sic, dtime, &
349         tsurf_cpl, p1lay, cal, beta, cdragh, cdragh, ps, &
350         precip_rain, precip_snow, snow, qsurf,  &
351         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
352         f_qsat_oce,AcoefH, AcoefQ, BcoefH, BcoefQ, &
353         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
354         sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol, rhoa)
355    do j = 1, knon
356      i = knindex(j)
357      sens_prec_liq_o(i,2) = sens_prec_liq(j)
358      sens_prec_sol_o(i,2) = sens_prec_sol(j)
359      lat_prec_liq_o(i,2) = lat_prec_liq(j)
360      lat_prec_sol_o(i,2) = lat_prec_sol(j)
361    enddo
362
363
364! - Flux calculation at first modele level for U and V
365    CALL calcul_flux_wind(knon, dtime, &
366         u0, v0, u1, v1, gustiness, cdragm, &
367         AcoefU, AcoefV, BcoefU, BcoefV, &
368         p1lay, temp_air, &
369         flux_u1, flux_v1) 
370
371!****************************************************************************************
372! Calculate fder : flux derivative (sensible and latente)
373!
374!****************************************************************************************
375    fder_new(:) = fder_old(:) + dflux_s(:) + dflux_l(:)
376   
377    iloc = MAXLOC(fder_new(1:klon))
378    IF (check .AND. fder_new(iloc(1))> 0.) THEN
379       WRITE(*,*)'**** Debug fder ****'
380       WRITE(*,*)'max fder(',iloc(1),') = ',fder_new(iloc(1))
381       WRITE(*,*)'fder_old, dflux_s, dflux_l',fder_old(iloc(1)), &
382            dflux_s(iloc(1)), dflux_l(iloc(1))
383    ENDIF
384
385!****************************************************************************************
386! Send and cumulate fields to the coupler
387!
388!****************************************************************************************
389
390    CALL cpl_send_seaice_fields(itime, dtime, knon, knindex, &
391       pctsrf, lafin, rlon, rlat, &
392       swnet, lwnet, fluxlat, fluxsens, &
393       precip_rain, precip_snow, evap, tsurf_new, fder_new, alb1, flux_u1, flux_v1,&
394       sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol)
395
396 
397
398  END SUBROUTINE ocean_cpl_ice
399
400!****************************************************************************************
401!
402END MODULE ocean_cpl_mod
Note: See TracBrowser for help on using the repository browser.