[781] | 1 | ! |
---|
[2538] | 2 | ! $Id: surf_ocean_mod.F90 2538 2016-06-03 14:12:16Z jyg $ |
---|
| 3 | ! |
---|
[781] | 4 | MODULE surf_ocean_mod |
---|
| 5 | |
---|
| 6 | IMPLICIT NONE |
---|
| 7 | |
---|
| 8 | CONTAINS |
---|
| 9 | ! |
---|
[2254] | 10 | !****************************************************************************** |
---|
[781] | 11 | ! |
---|
[888] | 12 | SUBROUTINE surf_ocean(rlon, rlat, swnet, lwnet, alb1, & |
---|
[2243] | 13 | windsp, rmu0, fder, tsurf_in, & |
---|
[781] | 14 | itime, dtime, jour, knon, knindex, & |
---|
[2254] | 15 | p1lay, z1lay, cdragh, cdragm, precip_rain, precip_snow, temp_air, spechum, & |
---|
[1067] | 16 | AcoefH, AcoefQ, BcoefH, BcoefQ, & |
---|
| 17 | AcoefU, AcoefV, BcoefU, BcoefV, & |
---|
[2240] | 18 | ps, u1, v1, gustiness, rugoro, pctsrf, & |
---|
[888] | 19 | snow, qsurf, agesno, & |
---|
[2243] | 20 | z0m, z0h, SFRWL, alb_dir_new, alb_dif_new, evap, fluxsens, fluxlat, & |
---|
[1067] | 21 | tsurf_new, dflux_s, dflux_l, lmt_bils, & |
---|
| 22 | flux_u1, flux_v1) |
---|
| 23 | |
---|
[2322] | 24 | use albedo, only: alboc, alboc_cd |
---|
[2391] | 25 | USE dimphy, ONLY: klon, zmasq |
---|
[1067] | 26 | USE surface_data, ONLY : type_ocean |
---|
| 27 | USE ocean_forced_mod, ONLY : ocean_forced_noice |
---|
| 28 | USE ocean_slab_mod, ONLY : ocean_slab_noice |
---|
| 29 | USE ocean_cpl_mod, ONLY : ocean_cpl_noice |
---|
[2391] | 30 | USE indice_sol_mod, ONLY : nbsrf, is_oce |
---|
[781] | 31 | ! |
---|
| 32 | ! This subroutine will make a call to ocean_XXX_noice according to the ocean mode (force, |
---|
| 33 | ! slab or couple). The calculations of albedo and rugosity for the ocean surface are |
---|
| 34 | ! done in here because they are identical for the different modes of ocean. |
---|
[1785] | 35 | |
---|
| 36 | |
---|
[793] | 37 | INCLUDE "YOMCST.h" |
---|
[781] | 38 | |
---|
[2178] | 39 | include "clesphys.h" |
---|
[2455] | 40 | ! for cycle_diurne and for iflag_z0_oce==-1 (prescribed z0) |
---|
[2178] | 41 | |
---|
[781] | 42 | ! Input variables |
---|
[2254] | 43 | !****************************************************************************** |
---|
[781] | 44 | INTEGER, INTENT(IN) :: itime, jour, knon |
---|
| 45 | INTEGER, DIMENSION(klon), INTENT(IN) :: knindex |
---|
| 46 | REAL, INTENT(IN) :: dtime |
---|
| 47 | REAL, DIMENSION(klon), INTENT(IN) :: rlon, rlat |
---|
[888] | 48 | REAL, DIMENSION(klon), INTENT(IN) :: swnet ! net shortwave radiation at surface |
---|
| 49 | REAL, DIMENSION(klon), INTENT(IN) :: lwnet ! net longwave radiation at surface |
---|
| 50 | REAL, DIMENSION(klon), INTENT(IN) :: alb1 ! albedo in visible SW interval |
---|
[781] | 51 | REAL, DIMENSION(klon), INTENT(IN) :: windsp |
---|
| 52 | REAL, DIMENSION(klon), INTENT(IN) :: rmu0 |
---|
| 53 | REAL, DIMENSION(klon), INTENT(IN) :: fder |
---|
[996] | 54 | REAL, DIMENSION(klon), INTENT(IN) :: tsurf_in |
---|
[2254] | 55 | REAL, DIMENSION(klon), INTENT(IN) :: p1lay,z1lay ! pression (Pa) et altitude (m) du premier niveau |
---|
[1067] | 56 | REAL, DIMENSION(klon), INTENT(IN) :: cdragh |
---|
| 57 | REAL, DIMENSION(klon), INTENT(IN) :: cdragm |
---|
[781] | 58 | REAL, DIMENSION(klon), INTENT(IN) :: precip_rain, precip_snow |
---|
| 59 | REAL, DIMENSION(klon), INTENT(IN) :: temp_air, spechum |
---|
[1067] | 60 | REAL, DIMENSION(klon), INTENT(IN) :: AcoefH, AcoefQ, BcoefH, BcoefQ |
---|
| 61 | REAL, DIMENSION(klon), INTENT(IN) :: AcoefU, AcoefV, BcoefU, BcoefV |
---|
[781] | 62 | REAL, DIMENSION(klon), INTENT(IN) :: ps |
---|
[2240] | 63 | REAL, DIMENSION(klon), INTENT(IN) :: u1, v1, gustiness |
---|
[781] | 64 | REAL, DIMENSION(klon), INTENT(IN) :: rugoro |
---|
| 65 | REAL, DIMENSION(klon,nbsrf), INTENT(IN) :: pctsrf |
---|
| 66 | |
---|
| 67 | ! In/Output variables |
---|
[2254] | 68 | !****************************************************************************** |
---|
[888] | 69 | REAL, DIMENSION(klon), INTENT(INOUT) :: snow |
---|
| 70 | REAL, DIMENSION(klon), INTENT(INOUT) :: qsurf |
---|
[781] | 71 | REAL, DIMENSION(klon), INTENT(INOUT) :: agesno |
---|
| 72 | |
---|
| 73 | ! Output variables |
---|
[2254] | 74 | !****************************************************************************** |
---|
[2243] | 75 | REAL, DIMENSION(klon), INTENT(OUT) :: z0m, z0h |
---|
[2227] | 76 | !albedo SB >>> |
---|
| 77 | ! REAL, DIMENSION(klon), INTENT(OUT) :: alb1_new ! new albedo in visible SW interval |
---|
| 78 | ! REAL, DIMENSION(klon), INTENT(OUT) :: alb2_new ! new albedo in near IR interval |
---|
| 79 | REAL, DIMENSION(6), INTENT(IN) :: SFRWL |
---|
| 80 | REAL, DIMENSION(klon,nsw), INTENT(OUT) :: alb_dir_new,alb_dif_new |
---|
| 81 | !albedo SB <<< |
---|
[781] | 82 | REAL, DIMENSION(klon), INTENT(OUT) :: evap, fluxsens, fluxlat |
---|
[888] | 83 | REAL, DIMENSION(klon), INTENT(OUT) :: tsurf_new |
---|
[781] | 84 | REAL, DIMENSION(klon), INTENT(OUT) :: dflux_s, dflux_l |
---|
[996] | 85 | REAL, DIMENSION(klon), INTENT(OUT) :: lmt_bils |
---|
[1067] | 86 | REAL, DIMENSION(klon), INTENT(OUT) :: flux_u1, flux_v1 |
---|
[781] | 87 | |
---|
| 88 | ! Local variables |
---|
[2254] | 89 | !****************************************************************************** |
---|
[2227] | 90 | INTEGER :: i, k |
---|
[1146] | 91 | REAL :: tmp |
---|
| 92 | REAL, PARAMETER :: cepdu2=(0.1)**2 |
---|
[781] | 93 | REAL, DIMENSION(klon) :: alb_eau |
---|
[888] | 94 | REAL, DIMENSION(klon) :: radsol |
---|
[2254] | 95 | REAL, DIMENSION(klon) :: cdragq ! Cdrag pour l'evaporation |
---|
[2391] | 96 | CHARACTER(len=20),PARAMETER :: modname="surf_ocean" |
---|
[781] | 97 | |
---|
| 98 | ! End definition |
---|
[2254] | 99 | !****************************************************************************** |
---|
[888] | 100 | |
---|
| 101 | |
---|
[2254] | 102 | !****************************************************************************** |
---|
[888] | 103 | ! Calculate total net radiance at surface |
---|
| 104 | ! |
---|
[2254] | 105 | !****************************************************************************** |
---|
[888] | 106 | radsol(:) = 0.0 |
---|
| 107 | radsol(1:knon) = swnet(1:knon) + lwnet(1:knon) |
---|
| 108 | |
---|
[2254] | 109 | !****************************************************************************** |
---|
| 110 | ! Cdragq computed from cdrag |
---|
| 111 | ! The difference comes only from a factor (f_z0qh_oce) on z0, so that |
---|
| 112 | ! it can be computed inside surf_ocean |
---|
| 113 | ! More complicated appraches may require the propagation through |
---|
| 114 | ! pbl_surface of an independant cdragq variable. |
---|
| 115 | !****************************************************************************** |
---|
| 116 | |
---|
| 117 | IF ( f_z0qh_oce .ne. 1.) THEN |
---|
[2261] | 118 | ! Si on suit les formulations par exemple de Tessel, on |
---|
| 119 | ! a z0h=0.4*nu/u*, z0q=0.62*nu/u*, d'ou f_z0qh_oce=0.62/0.4=1.55 |
---|
[2254] | 120 | cdragq(:)=cdragh(:)* & |
---|
| 121 | log(z1lay(:)/z0h(:))/log(z1lay(:)/(f_z0qh_oce*z0h(:))) |
---|
| 122 | ELSE |
---|
| 123 | cdragq(:)=cdragh(:) |
---|
| 124 | ENDIF |
---|
| 125 | |
---|
| 126 | !****************************************************************************** |
---|
[781] | 127 | ! Switch according to type of ocean (couple, slab or forced) |
---|
[2254] | 128 | !****************************************************************************** |
---|
[996] | 129 | SELECT CASE(type_ocean) |
---|
[781] | 130 | CASE('couple') |
---|
[888] | 131 | CALL ocean_cpl_noice( & |
---|
| 132 | swnet, lwnet, alb1, & |
---|
[1067] | 133 | windsp, fder, & |
---|
[781] | 134 | itime, dtime, knon, knindex, & |
---|
[2254] | 135 | p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow,temp_air,spechum,& |
---|
[1067] | 136 | AcoefH, AcoefQ, BcoefH, BcoefQ, & |
---|
| 137 | AcoefU, AcoefV, BcoefU, BcoefV, & |
---|
[2240] | 138 | ps, u1, v1, gustiness, & |
---|
[888] | 139 | radsol, snow, agesno, & |
---|
[1067] | 140 | qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, & |
---|
[996] | 141 | tsurf_new, dflux_s, dflux_l) |
---|
[781] | 142 | |
---|
| 143 | CASE('slab') |
---|
[888] | 144 | CALL ocean_slab_noice( & |
---|
[996] | 145 | itime, dtime, jour, knon, knindex, & |
---|
[2254] | 146 | p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, temp_air, spechum,& |
---|
[1067] | 147 | AcoefH, AcoefQ, BcoefH, BcoefQ, & |
---|
| 148 | AcoefU, AcoefV, BcoefU, BcoefV, & |
---|
[2240] | 149 | ps, u1, v1, gustiness, tsurf_in, & |
---|
[2209] | 150 | radsol, snow, & |
---|
[1067] | 151 | qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, & |
---|
[996] | 152 | tsurf_new, dflux_s, dflux_l, lmt_bils) |
---|
[781] | 153 | |
---|
| 154 | CASE('force') |
---|
[888] | 155 | CALL ocean_forced_noice( & |
---|
| 156 | itime, dtime, jour, knon, knindex, & |
---|
[2254] | 157 | p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, & |
---|
[781] | 158 | temp_air, spechum, & |
---|
[1067] | 159 | AcoefH, AcoefQ, BcoefH, BcoefQ, & |
---|
| 160 | AcoefU, AcoefV, BcoefU, BcoefV, & |
---|
[2240] | 161 | ps, u1, v1, gustiness, & |
---|
[888] | 162 | radsol, snow, agesno, & |
---|
[1067] | 163 | qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, & |
---|
[996] | 164 | tsurf_new, dflux_s, dflux_l) |
---|
[781] | 165 | END SELECT |
---|
| 166 | |
---|
[2254] | 167 | !****************************************************************************** |
---|
[2057] | 168 | ! fcodron: compute lmt_bils forced case (same as wfbils_oce / 1.-contfracatm) |
---|
[2254] | 169 | !****************************************************************************** |
---|
[2057] | 170 | IF (type_ocean.NE.'slab') THEN |
---|
| 171 | lmt_bils(:)=0. |
---|
| 172 | DO i=1,knon |
---|
| 173 | lmt_bils(knindex(i))=(swnet(i)+lwnet(i)+fluxsens(i)+fluxlat(i)) & |
---|
| 174 | *pctsrf(knindex(i),is_oce)/(1.-zmasq(knindex(i))) |
---|
| 175 | END DO |
---|
| 176 | END IF |
---|
| 177 | |
---|
[2254] | 178 | !****************************************************************************** |
---|
[2413] | 179 | ! Calculate ocean surface albedo |
---|
[2254] | 180 | !****************************************************************************** |
---|
[2227] | 181 | !albedo SB >>> |
---|
[2413] | 182 | IF (iflag_albedo==0) THEN |
---|
| 183 | !--old parametrizations of ocean surface albedo |
---|
| 184 | ! |
---|
[2178] | 185 | IF (cycle_diurne) THEN |
---|
[2413] | 186 | ! |
---|
[2178] | 187 | CALL alboc_cd(rmu0,alb_eau) |
---|
[2413] | 188 | ! |
---|
| 189 | !--ad-hoc correction for model radiative balance tuning |
---|
| 190 | !--now outside alboc_cd routine |
---|
| 191 | alb_eau(:) = fmagic*alb_eau(:) + pmagic |
---|
| 192 | alb_eau=MIN(MAX(alb_eau,0.0),1.0) |
---|
| 193 | ! |
---|
[2178] | 194 | ELSE |
---|
[2413] | 195 | ! |
---|
[1403] | 196 | CALL alboc(REAL(jour),rlat,alb_eau) |
---|
[2413] | 197 | !--ad-hoc correction for model radiative balance tuning |
---|
| 198 | !--now outside alboc routine |
---|
| 199 | alb_eau(:) = fmagic*alb_eau(:) + pmagic |
---|
| 200 | alb_eau=MIN(MAX(alb_eau(i),0.04),0.60) |
---|
| 201 | ! |
---|
[781] | 202 | ENDIF |
---|
[2413] | 203 | ! |
---|
[781] | 204 | DO i =1, knon |
---|
[2413] | 205 | DO k=1,nsw |
---|
[2227] | 206 | alb_dir_new(i,k) = alb_eau(knindex(i)) |
---|
[2413] | 207 | ENDDO |
---|
[781] | 208 | ENDDO |
---|
[2413] | 209 | !IM 09122015 next line corresponds to the old way of doing in LMDZ5A/IPSLCM5A versions |
---|
| 210 | !albedo for diffuse radiation is taken the same as for direct radiation |
---|
[2405] | 211 | alb_dif_new=alb_dir_new |
---|
| 212 | !IM 09122015 end |
---|
[2413] | 213 | ! |
---|
| 214 | ELSE IF (iflag_albedo==1) THEN |
---|
| 215 | !--new parametrization of ocean surface albedo by Sunghye Baek |
---|
| 216 | !--albedo for direct and diffuse radiation are different |
---|
| 217 | ! |
---|
| 218 | CALL ocean_albedo(knon,rmu0,knindex,windsp,SFRWL,alb_dir_new,alb_dif_new) |
---|
| 219 | ! |
---|
| 220 | !--ad-hoc correction for model radiative balance tuning |
---|
| 221 | alb_dir_new(:,:) = fmagic*alb_dir_new(:,:) + pmagic |
---|
| 222 | alb_dir_new=MIN(MAX(alb_dir_new,0.0),1.0) |
---|
| 223 | alb_dif_new=MIN(MAX(alb_dif_new,0.0),1.0) |
---|
| 224 | ! |
---|
| 225 | ENDIF |
---|
[2227] | 226 | !albedo SB <<< |
---|
| 227 | |
---|
[2254] | 228 | !****************************************************************************** |
---|
[781] | 229 | ! Calculate the rugosity |
---|
[2254] | 230 | !****************************************************************************** |
---|
[2243] | 231 | IF (iflag_z0_oce==0) THEN |
---|
[781] | 232 | DO i = 1, knon |
---|
[2278] | 233 | tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2) |
---|
| 234 | z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG & |
---|
[1146] | 235 | + 0.11*14e-6 / SQRT(cdragm(i) * tmp) |
---|
[2243] | 236 | z0m(i) = MAX(1.5e-05,z0m(i)) |
---|
[1146] | 237 | ENDDO |
---|
[2243] | 238 | z0h(1:knon)=z0m(1:knon) ! En attendant mieux |
---|
| 239 | |
---|
[2261] | 240 | ELSE IF (iflag_z0_oce==1) THEN |
---|
[2264] | 241 | DO i = 1, knon |
---|
[2278] | 242 | tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2) |
---|
| 243 | z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG & |
---|
[2264] | 244 | + 0.11*14e-6 / SQRT(cdragm(i) * tmp) |
---|
| 245 | z0m(i) = MAX(1.5e-05,z0m(i)) |
---|
[2278] | 246 | z0h(i)=0.4*14e-6 / SQRT(cdragm(i) * tmp) |
---|
[2264] | 247 | ENDDO |
---|
[2455] | 248 | ELSE IF (iflag_z0_oce==-1) THEN |
---|
| 249 | DO i = 1, knon |
---|
| 250 | z0m(i) = z0min |
---|
| 251 | z0h(i) = z0min |
---|
| 252 | ENDDO |
---|
[2243] | 253 | ELSE |
---|
[2391] | 254 | CALL abort_physic(modname,'version non prevue',1) |
---|
[2243] | 255 | ENDIF |
---|
[781] | 256 | ! |
---|
[2254] | 257 | !****************************************************************************** |
---|
[781] | 258 | END SUBROUTINE surf_ocean |
---|
[2254] | 259 | !****************************************************************************** |
---|
[781] | 260 | ! |
---|
| 261 | END MODULE surf_ocean_mod |
---|