Changeset 1761 for LMDZ5/trunk
- Timestamp:
- Jun 4, 2013, 12:11:38 PM (11 years ago)
- Location:
- LMDZ5/trunk/libf/phylmd
- Files:
-
- 1 added
- 10 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/phylmd/coef_diff_turb_mod.F90
r1738 r1761 158 158 159 159 ! iflag_pbl peut etre utilise comme longuer de melange 160 IF (iflag_pbl.GE. 18) THEN160 IF (iflag_pbl.GE.31) THEN 161 161 CALL vdif_kcay(knon,dtime,RG,RD,ypaprs,yt, & 162 162 yzlev,yzlay,yu,yv,yteta, & 163 163 ycdragm,yq2,q2diag,ykmm,ykmn,yustar, & 164 164 iflag_pbl) 165 ELSE 165 ELSE IF (iflag_pbl<20) THEN 166 166 CALL yamada4(knon,dtime,RG,RD,ypaprs,yt, & 167 167 yzlev,yzlay,yu,yv,yteta, & -
LMDZ5/trunk/libf/phylmd/ener_conserv.F90
r1754 r1761 1 1 subroutine ener_conserv(klon,klev,pdtphys, & 2 & puo,pvo,pto,pqo,pun,pvn,ptn,pqn, masse,exner,d_t_ec)2 & puo,pvo,pto,pqo,pun,pvn,ptn,pqn,dtke,masse,exner,d_t_ec) 3 3 4 4 !============================================================= … … 19 19 20 20 ! From module 21 USE phys_local_var_mod, ONLY : d_u_vdf,d_v_vdf,d_t_vdf,d_u_ajs,d_v_ajs,d_t_ajs,d_u_con,d_v_con,d_t_con 22 USE phys_output_var_mod, ONLY : bils_ec,bils_ kinetic,bils_enthalp,bils_latent21 USE phys_local_var_mod, ONLY : d_u_vdf,d_v_vdf,d_t_vdf,d_u_ajs,d_v_ajs,d_t_ajs,d_u_con,d_v_con,d_t_con,d_t_diss 22 USE phys_output_var_mod, ONLY : bils_ec,bils_tke,bils_kinetic,bils_enthalp,bils_latent,bils_diss 23 23 24 24 IMPLICIT none … … 26 26 #include "YOETHF.h" 27 27 #include "clesphys.h" 28 #include "compbl.h" 28 29 29 30 ! Arguments … … 33 34 REAL, DIMENSION(klon,klev),INTENT(IN) :: pun,pvn,ptn,pqn 34 35 REAL, DIMENSION(klon,klev),INTENT(IN) :: masse,exner 36 REAL, DIMENSION(klon,klev+1),INTENT(IN) :: dtke 35 37 REAL, DIMENSION(klon,klev),INTENT(OUT) :: d_t_ec 36 38 integer k,i … … 65 67 66 68 IF (iflag_ener_conserv<=2) THEN 67 d_t(:,:)=d_t_vdf(:,:)+d_t_ajs(:,:) ! d_t_ajs = adjust + thermals 68 d_u(:,:)=d_u_vdf(:,:)+d_u_ajs(:,:)+d_u_con(:,:) 69 d_v(:,:)=d_v_vdf(:,:)+d_v_ajs(:,:)+d_v_con(:,:) 69 ! print*,'ener_conserv pbl=',iflag_pbl 70 IF (iflag_pbl>=20 .AND. iflag_pbl<=27) THEN !d_t_diss accounts for conserv 71 d_t(:,:)=d_t_ajs(:,:) ! d_t_ajs = adjust + thermals 72 d_u(:,:)=d_u_ajs(:,:)+d_u_con(:,:) 73 d_v(:,:)=d_v_ajs(:,:)+d_v_con(:,:) 74 ELSE 75 d_t(:,:)=d_t_vdf(:,:)+d_t_ajs(:,:) ! d_t_ajs = adjust + thermals 76 d_u(:,:)=d_u_vdf(:,:)+d_u_ajs(:,:)+d_u_con(:,:) 77 d_v(:,:)=d_v_vdf(:,:)+d_v_ajs(:,:)+d_v_con(:,:) 78 ENDIF 70 79 ELSEIF (iflag_ener_conserv==101) THEN 71 80 d_t(:,:)=0. … … 89 98 zv(:,:)=pvo(:,:) 90 99 else 91 zu(:,:)=puo(:,:)+0.5*d_u(:,:) 92 zv(:,:)=pvo(:,:)+0.5*d_v(:,:) 100 IF (iflag_pbl>=20 .AND. iflag_pbl<=27) THEN 101 zu(:,:)=puo(:,:)+d_u_vdf(:,:)+0.5*d_u(:,:) 102 zv(:,:)=pvo(:,:)+d_v_vdf(:,:)+0.5*d_v(:,:) 103 ELSE 104 zu(:,:)=puo(:,:)+0.5*d_u(:,:) 105 zv(:,:)=pvo(:,:)+0.5*d_v(:,:) 106 ENDIF 93 107 endif 94 108 … … 130 144 131 145 bils_ec(:)=0. 146 bils_tke(:)=0. 147 bils_diss(:)=0. 132 148 bils_kinetic(:)=0. 133 149 bils_enthalp(:)=0. … … 135 151 DO k=1,klev 136 152 bils_ec(:)=bils_ec(:)-d_t_ec(:,k)*masse(:,k) 153 bils_tke(:)=bils_tke(:)+0.5*(dtke(:,k)+dtke(:,k+1))*masse(:,k) 154 bils_diss(:)=bils_diss(:)-d_t_diss(:,k)*masse(:,k) 137 155 bils_kinetic(:)=bils_kinetic(:)+masse(:,k)* & 138 156 & (pun(:,k)*pun(:,k)+pvn(:,k)*pvn(:,k) & … … 144 162 ENDDO 145 163 bils_ec(:)=rcpd*bils_ec(:)/pdtphys 164 bils_tke(:)=bils_tke(:)/pdtphys 165 bils_diss(:)=rcpd*bils_diss(:)/pdtphys 146 166 bils_kinetic(:)= 0.5*bils_kinetic(:)/pdtphys 147 167 bils_enthalp(:)=rcpd*bils_enthalp(:)/pdtphys 148 168 bils_latent(:)=rlvtt*bils_latent(:)/pdtphys 149 150 151 169 RETURN 152 170 -
LMDZ5/trunk/libf/phylmd/indicesol.h
r1279 r1761 21 21 PARAMETER (is_lic=2) ! glacier continental 22 22 ! 23 INTEGER is_ave 24 PARAMETER (is_ave=nbsrf+1) ! glacier continental 25 ! 23 26 REAL epsfra 24 27 PARAMETER (epsfra=1.0E-05) -
LMDZ5/trunk/libf/phylmd/pbl_surface_mod.F90
r1670 r1761 176 176 alb1_m, alb2_m, zxsens, zxevap, & 177 177 zxtsol, zxfluxlat, zt2m, qsat2m, & 178 d_t, d_q, d_u, d_v, 178 d_t, d_q, d_u, d_v, d_t_diss, & 179 179 zcoefh, zcoefm, slab_wfbils, & 180 180 qsol_d, zq2m, s_pblh, s_plcl, & … … 259 259 INCLUDE "YOETHF.h" 260 260 INCLUDE "temps.h" 261 !**************************************************************************************** 262 ! Declarations specifiques pour le 1D. A reprendre 261 263 ! Input variables 262 264 !**************************************************************************************** … … 291 293 REAL, DIMENSION(klon, nbsrf), INTENT(INOUT) :: u10m ! u speed at 10m 292 294 REAL, DIMENSION(klon, nbsrf), INTENT(INOUT) :: v10m ! v speed at 10m 293 REAL, DIMENSION(klon, klev+1, nbsrf ), INTENT(INOUT) :: tke295 REAL, DIMENSION(klon, klev+1, nbsrf+1), INTENT(INOUT) :: tke 294 296 295 297 ! Output variables … … 310 312 REAL, DIMENSION(klon), INTENT(OUT) :: qsat2m 311 313 REAL, DIMENSION(klon, klev), INTENT(OUT) :: d_t ! change in temperature 314 REAL, DIMENSION(klon, klev), INTENT(OUT) :: d_t_diss ! change in temperature 312 315 REAL, DIMENSION(klon, klev), INTENT(OUT) :: d_q ! change in water vapour 313 316 REAL, DIMENSION(klon, klev), INTENT(OUT) :: d_u ! change in u speed 314 317 REAL, DIMENSION(klon, klev), INTENT(OUT) :: d_v ! change in v speed 315 REAL, DIMENSION(klon, klev ), INTENT(OUT) :: zcoefh ! coef for turbulent diffusion of T and Q, mean for each grid point316 REAL, DIMENSION(klon, klev ), INTENT(OUT) :: zcoefm ! coef for turbulent diffusion of U and V (?), mean for each grid point318 REAL, DIMENSION(klon, klev,nbsrf+1), INTENT(OUT) :: zcoefh ! coef for turbulent diffusion of T and Q, mean for each grid point 319 REAL, DIMENSION(klon, klev,nbsrf+1), INTENT(OUT) :: zcoefm ! coef for turbulent diffusion of U and V (?), mean for each grid point 317 320 318 321 ! Output only for diagnostics … … 424 427 REAL, DIMENSION(klon) :: ztsol 425 428 REAL, DIMENSION(klon) :: alb_m ! mean albedo for whole SW interval 426 REAL, DIMENSION(klon,klev) :: y_d_t, y_d_q 429 REAL, DIMENSION(klon,klev) :: y_d_t, y_d_q, y_d_t_diss 427 430 REAL, DIMENSION(klon,klev) :: y_d_u, y_d_v 428 431 REAL, DIMENSION(klon,klev) :: y_flux_t, y_flux_q 429 432 REAL, DIMENSION(klon,klev) :: y_flux_u, y_flux_v 430 REAL, DIMENSION(klon,klev) :: ycoefh, ycoefm 433 REAL, DIMENSION(klon,klev) :: ycoefh, ycoefm,ycoefq 431 434 REAL, DIMENSION(klon) :: ycdragh, ycdragm 432 435 REAL, DIMENSION(klon,klev) :: yu, yv … … 470 473 !**************************************************************************************** 471 474 ! Declarations specifiques pour le 1D. A reprendre 475 !**************************************************************************************** 472 476 REAL :: fsens,flat 473 477 LOGICAL :: ok_flux_surf ! initialized during first_call below 474 478 COMMON /flux_arp/fsens,flat,ok_flux_surf 475 !****************************************************************************************476 479 ! End of declarations 477 480 !**************************************************************************************** … … 546 549 d_ts = 0.0 ; yfluxlat=0.0 ; flux_t = 0.0 ; flux_q = 0.0 547 550 flux_u = 0.0 ; flux_v = 0.0 ; d_t = 0.0 ; d_q = 0.0 548 d_ u = 0.0 ; d_v = 0.0 ; yqsol = 0.0551 d_t_diss= 0.0 ;d_u = 0.0 ; d_v = 0.0 ; yqsol = 0.0 549 552 ytherm = 0.0 ; ytke=0. 550 553 551 zcoefh(:,:) = 0.0 552 zcoefh(:,1) = 999999. ! zcoefh(:,k=1) should never be used 553 zcoefm(:,:) = 0.0 554 zcoefm(:,1) = 999999. ! 554 tke(:,:,is_ave)=0. 555 IF (iflag_pbl<20.or.iflag_pbl>=30) THEN 556 zcoefh(:,:,:) = 0.0 557 zcoefh(:,1,:) = 999999. ! zcoefh(:,k=1) should never be used 558 zcoefm(:,:,:) = 0.0 559 zcoefm(:,1,:) = 999999. ! 560 ELSE 561 zcoefm(:,:,is_ave)=0. 562 zcoefh(:,:,is_ave)=0. 563 ENDIF 555 564 ytsoil = 999999. 556 565 … … 713 722 ENDDO 714 723 ENDDO 715 724 716 725 DO k = 1, nsoilmx 717 726 DO j = 1, knon … … 747 756 ypaprs, ypplay, yu, yv, yq, yt, yts, yrugos, yqsurf, ycdragm, & 748 757 ycoefm, ycoefh, ytke) 758 759 IF (iflag_pbl>=20.AND.iflag_pbl<30) THEN 760 ! In this case, coef_diff_turb is called for the Cd only 761 DO k = 2, klev 762 DO j = 1, knon 763 i = ni(j) 764 ycoefh(j,k) = zcoefh(i,k,nsrf) 765 ycoefm(j,k) = zcoefm(i,k,nsrf) 766 ENDDO 767 ENDDO 768 ENDIF 749 769 750 770 !**************************************************************************************** … … 924 944 925 945 946 y_d_t_diss(:,:)=0. 947 IF (iflag_pbl>=20 .and. iflag_pbl<30) THEN 948 CALL yamada_c(knon,dtime,ypaprs,ypplay & 949 & ,yu,yv,yt,y_d_u,y_d_v,y_d_t,ycdragm,ytke,ycoefm,ycoefh,ycoefq,y_d_t_diss,yustar & 950 & ,iflag_pbl,nsrf) 951 ENDIF 952 ! print*,'yamada_c OK' 953 926 954 DO j = 1, knon 927 955 y_dflux_t(j) = y_dflux_t(j) * ypct(j) … … 937 965 !**************************************************************************************** 938 966 939 tke(:,:,nsrf) = 0.940 967 DO k = 1, klev 941 968 DO j = 1, knon 942 969 i = ni(j) 970 y_d_t_diss(j,k) = y_d_t_diss(j,k) * ypct(j) 943 971 y_d_t(j,k) = y_d_t(j,k) * ypct(j) 944 972 y_d_q(j,k) = y_d_q(j,k) * ypct(j) … … 951 979 flux_v(i,k,nsrf) = y_flux_v(j,k) 952 980 953 tke(i,k,nsrf) = ytke(j,k)954 981 955 982 ENDDO 956 983 ENDDO 984 985 ! print*,'Dans pbl OK1' 957 986 958 987 evap(:,nsrf) = - flux_q(:,1,nsrf) … … 980 1009 END DO 981 1010 1011 ! print*,'Dans pbl OK2' 1012 982 1013 DO k = 2, klev 983 1014 DO j = 1, knon 984 1015 i = ni(j) 985 zcoefh(i,k) = zcoefh(i,k) + ycoefh(j,k)*ypct(j) 986 zcoefm(i,k) = zcoefm(i,k) + ycoefm(j,k)*ypct(j) 1016 tke(i,k,nsrf) = ytke(j,k) 1017 zcoefh(i,k,nsrf) = ycoefh(j,k) 1018 zcoefm(i,k,nsrf) = ycoefm(j,k) 1019 tke(i,k,is_ave) = tke(i,k,is_ave) + ytke(j,k)*ypct(j) 1020 zcoefh(i,k,is_ave) = zcoefh(i,k,is_ave) + ycoefh(j,k)*ypct(j) 1021 zcoefm(i,k,is_ave) = zcoefm(i,k,is_ave) + ycoefm(j,k)*ypct(j) 987 1022 END DO 988 1023 END DO 1024 1025 ! print*,'Dans pbl OK3' 989 1026 990 1027 IF ( nsrf .EQ. is_ter ) THEN … … 1007 1044 DO j = 1, knon 1008 1045 i = ni(j) 1046 d_t_diss(i,k) = d_t_diss(i,k) + y_d_t_diss(j,k) 1009 1047 d_t(i,k) = d_t(i,k) + y_d_t(j,k) 1010 1048 d_q(i,k) = d_q(i,k) + y_d_q(j,k) … … 1013 1051 END DO 1014 1052 END DO 1053 1054 ! print*,'Dans pbl OK4' 1015 1055 1016 1056 !**************************************************************************************** … … 1040 1080 ! Calculations of diagnostic t,q at 2m and u, v at 10m 1041 1081 1082 ! print*,'Dans pbl OK41' 1083 ! print*,'tair1,yt(:,1),y_d_t(:,1)' 1084 ! print*, tair1,yt(:,1),y_d_t(:,1) 1042 1085 DO j=1, knon 1043 1086 i = ni(j) 1044 1087 uzon(j) = yu(j,1) + y_d_u(j,1) 1045 1088 vmer(j) = yv(j,1) + y_d_v(j,1) 1046 tair1(j) = yt(j,1) + y_d_t(j,1) 1089 tair1(j) = yt(j,1) + y_d_t(j,1) + y_d_t_diss(j,1) 1047 1090 qair1(j) = yq(j,1) + y_d_q(j,1) 1048 1091 zgeo1(j) = RD * tair1(j) / (0.5*(ypaprs(j,1)+ypplay(j,1))) & … … 1058 1101 END DO 1059 1102 1103 ! print*,'Dans pbl OK42A' 1104 ! print*,'tair1,yt(:,1),y_d_t(:,1)' 1105 ! print*, tair1,yt(:,1),y_d_t(:,1) 1060 1106 1061 1107 ! Calculate the temperature et relative humidity at 2m and the wind at 10m … … 1064 1110 tairsol, qairsol, rugo1, psfce, patm, & 1065 1111 yt2m, yq2m, yt10m, yq10m, yu10m, yustar) 1112 ! print*,'Dans pbl OK42B' 1066 1113 1067 1114 DO j=1, knon … … 1077 1124 END DO 1078 1125 1126 ! print*,'Dans pbl OK43' 1079 1127 !IM Calcule de l'humidite relative a 2m (rh2m) pour diagnostique 1080 1128 !IM Ajoute dependance type surface … … 1093 1141 END IF 1094 1142 1143 ! print*,'OK pbl 5' 1095 1144 CALL HBTM(knon, ypaprs, ypplay, & 1096 1145 yt2m,yt10m,yq2m,yq10m,yustar, & … … 1113 1162 END DO 1114 1163 1164 ! print*,'OK pbl 6' 1115 1165 #else 1116 1166 ! T2m not defined … … 1130 1180 !**************************************************************************************** 1131 1181 1182 ! print*,'OK pbl 7' 1132 1183 zxfluxt(:,:) = 0.0 ; zxfluxq(:,:) = 0.0 1133 1184 zxfluxu(:,:) = 0.0 ; zxfluxv(:,:) = 0.0 … … 1143 1194 END DO 1144 1195 1196 ! print*,'OK pbl 8' 1145 1197 DO i = 1, klon 1146 1198 zxsens(i) = - zxfluxt(i,1) ! flux de chaleur sensible au sol … … 1161 1213 s_trmb2(:) = 0.0 ; s_trmb3(:) = 0.0 1162 1214 1215 ! print*,'OK pbl 9' 1163 1216 1164 1217 DO nsrf = 1, nbsrf … … 1192 1245 END DO 1193 1246 END DO 1247 ! print*,'OK pbl 10' 1194 1248 1195 1249 IF (check) THEN -
LMDZ5/trunk/libf/phylmd/phys_local_var_mod.F90
r1742 r1761 55 55 REAL, SAVE, ALLOCATABLE :: d_u_oli(:,:), d_v_oli(:,:) 56 56 !$OMP THREADPRIVATE(d_u_oli, d_v_oli) 57 REAL, SAVE, ALLOCATABLE :: d_t_vdf(:,:), d_q_vdf(:,:) 58 !$OMP THREADPRIVATE( d_t_vdf, d_q_vdf )57 REAL, SAVE, ALLOCATABLE :: d_t_vdf(:,:), d_q_vdf(:,:), d_t_diss(:,:) 58 !$OMP THREADPRIVATE( d_t_vdf, d_q_vdf,d_t_diss) 59 59 REAL, SAVE, ALLOCATABLE :: d_u_vdf(:,:), d_v_vdf(:,:) 60 60 !$OMP THREADPRIVATE(d_u_vdf, d_v_vdf) … … 216 216 allocate(d_t_lscth(klon,klev),d_q_lscth(klon,klev)) 217 217 allocate(plul_st(klon),plul_th(klon)) 218 allocate(d_t_vdf(klon,klev),d_q_vdf(klon,klev) )218 allocate(d_t_vdf(klon,klev),d_q_vdf(klon,klev),d_t_diss(klon,klev)) 219 219 allocate(d_u_vdf(klon,klev),d_v_vdf(klon,klev)) 220 220 allocate(d_t_oli(klon,klev),d_t_oro(klon,klev)) … … 305 305 deallocate(d_t_lscth,d_q_lscth) 306 306 deallocate(plul_st,plul_th) 307 deallocate(d_t_vdf,d_q_vdf )307 deallocate(d_t_vdf,d_q_vdf,d_t_diss) 308 308 deallocate(d_u_vdf,d_v_vdf) 309 309 deallocate(d_t_oli,d_t_oro) -
LMDZ5/trunk/libf/phylmd/phys_output_mod.F90
r1753 r1761 156 156 type(ctrl_out),save :: o_LWdnSFCclr = ctrl_out((/ 1, 4, 10, 10, 5, 10 /),'LWdnSFCclr') 157 157 type(ctrl_out),save :: o_bils = ctrl_out((/ 1, 2, 10, 5, 10, 10 /),'bils') 158 type(ctrl_out),save :: o_bils_tke = ctrl_out((/ 1, 2, 10, 5, 10, 10 /),'bils_tke') 159 type(ctrl_out),save :: o_bils_diss = ctrl_out((/ 1, 2, 10, 5, 10, 10 /),'bils_diss') 158 160 type(ctrl_out),save :: o_bils_ec = ctrl_out((/ 1, 2, 10, 5, 10, 10 /),'bils_ec') 159 161 type(ctrl_out),save :: o_bils_kinetic = ctrl_out((/ 1, 2, 10, 5, 10, 10 /),'bils_kinetic') … … 554 556 type(ctrl_out),save :: o_beta_prec = ctrl_out((/ 4, 10, 10, 10, 10, 10 /),'beta_prec') 555 557 type(ctrl_out),save :: o_dtvdf = ctrl_out((/ 4, 10, 10, 10, 10, 10 /),'dtvdf') 558 type(ctrl_out),save :: o_dtdis = ctrl_out((/ 4, 10, 10, 10, 10, 10 /),'dtdis') 556 559 type(ctrl_out),save :: o_dqvdf = ctrl_out((/ 4, 10, 10, 10, 10, 10 /),'dqvdf') 557 560 type(ctrl_out),save :: o_dteva = ctrl_out((/ 4, 10, 10, 10, 10, 10 /),'dteva') … … 1043 1046 CALL histdef2d(iff,clef_stations(iff),o_bils%flag,o_bils%name, "Surf. total heat flux", "W/m2") 1044 1047 CALL histdef2d(iff,clef_stations(iff),o_bils_ec%flag,o_bils_ec%name, "Surf. total heat flux", "W/m2") 1048 CALL histdef2d(iff,clef_stations(iff),o_bils_tke%flag,o_bils_tke%name, "Surf. total heat flux", "W/m2") 1049 CALL histdef2d(iff,clef_stations(iff),o_bils_diss%flag,o_bils_diss%name, "Surf. total heat flux", "W/m2") 1045 1050 CALL histdef2d(iff,clef_stations(iff),o_bils_kinetic%flag,o_bils_kinetic%name, "Surf. total heat flux", "W/m2") 1046 1051 CALL histdef2d(iff,clef_stations(iff),o_bils_enthalp%flag,o_bils_enthalp%name, "Surf. total heat flux", "W/m2") … … 1625 1630 CALL histdef3d(iff,clef_stations(iff),o_beta_prec%flag,o_beta_prec%name, "LS Conversion rate to prec", "(kg/kg)/s") 1626 1631 CALL histdef3d(iff,clef_stations(iff),o_dtvdf%flag,o_dtvdf%name, "Boundary-layer dT", "K/s") 1632 CALL histdef3d(iff,clef_stations(iff),o_dtdis%flag,o_dtdis%name, "TKE dissipation dT", "K/s") 1627 1633 CALL histdef3d(iff,clef_stations(iff),o_dqvdf%flag,o_dqvdf%name, "Boundary-layer dQ", "(kg/kg)/s") 1628 1634 CALL histdef3d(iff,clef_stations(iff),o_dteva%flag,o_dteva%name, "Reevaporation dT", "K/s") -
LMDZ5/trunk/libf/phylmd/phys_output_var_mod.F90
r1758 r1761 17 17 !$OMP THREADPRIVATE(itau_con) 18 18 REAL, ALLOCATABLE :: bils_ec(:) ! Contribution of energy conservation 19 REAL, ALLOCATABLE :: bils_tke(:) ! Contribution of energy conservation 20 REAL, ALLOCATABLE :: bils_diss(:) ! Contribution of energy conservation 19 21 REAL, ALLOCATABLE :: bils_kinetic(:) ! bilan de chaleur au sol, kinetic 20 22 REAL, ALLOCATABLE :: bils_enthalp(:) ! bilan de chaleur au sol 21 23 REAL, ALLOCATABLE :: bils_latent(:) ! bilan de chaleur au sol 22 !$OMP THREADPRIVATE(bils_ec,bils_kinetic,bils_enthalp,bils_latent) 24 !$OMP THREADPRIVATE(bils_ec,bils_tke,bils_diss,bils_kinetic,bils_enthalp,bils_latent) 25 23 26 24 27 CONTAINS … … 32 35 allocate(snow_o(klon), zfra_o(klon)) 33 36 allocate(itau_con(klon)) 34 allocate (bils_ec(klon),bils_ kinetic(klon),bils_enthalp(klon),bils_latent(klon))37 allocate (bils_ec(klon),bils_tke(klon),bils_diss(klon),bils_kinetic(klon),bils_enthalp(klon),bils_latent(klon)) 35 38 36 39 END SUBROUTINE phys_output_var_init … … 42 45 43 46 deallocate(snow_o,zfra_o,itau_con) 44 deallocate (bils_ec,bils_ kinetic,bils_enthalp,bils_latent)47 deallocate (bils_ec,bils_tke,bils_diss,bils_kinetic,bils_enthalp,bils_latent) 45 48 46 49 END SUBROUTINE phys_output_var_end -
LMDZ5/trunk/libf/phylmd/phys_output_write.h
r1753 r1761 360 360 ENDIF 361 361 362 IF (o_bils_diss%flag(iff)<=lev_files(iff)) THEN 363 CALL histwrite_phy(nid_files(iff),clef_stations(iff), 364 $o_bils_diss%name,itau_w,bils_diss) 365 ENDIF 366 362 367 IF (o_bils_ec%flag(iff)<=lev_files(iff)) THEN 363 368 CALL histwrite_phy(nid_files(iff),clef_stations(iff), 364 369 $o_bils_ec%name,itau_w,bils_ec) 370 ENDIF 371 372 IF (o_bils_tke%flag(iff)<=lev_files(iff)) THEN 373 CALL histwrite_phy(nid_files(iff),clef_stations(iff), 374 $o_bils_tke%name,itau_w,bils_tke) 365 375 ENDIF 366 376 … … 1333 1343 1334 1344 !====MS forcing diagnostics 1335 if (new_aod) then 1345 if (new_aod) then 1336 1346 IF (o_swtoaas_nat%flag(iff)<=lev_files(iff)) THEN 1337 1347 CALL histwrite_phy(nid_files(iff),clef_stations(iff), … … 1715 1725 IF (o_kz%flag(iff)<=lev_files(iff)) THEN 1716 1726 CALL histwrite_phy(nid_files(iff),clef_stations(iff), 1717 $o_kz%name,itau_w,coefh )1727 $o_kz%name,itau_w,coefh(:,:,is_ave)) 1718 1728 ENDIF 1719 1729 … … 1721 1731 IF (o_kz_max%flag(iff)<=lev_files(iff)) THEN 1722 1732 CALL histwrite_phy(nid_files(iff),clef_stations(iff), 1723 $o_kz_max%name,itau_w,coefh )1733 $o_kz_max%name,itau_w,coefh(:,:,is_ave)) 1724 1734 ENDIF 1725 1735 ENDIF … … 1908 1918 ENDIF 1909 1919 1920 IF (o_dtdis%flag(iff)<=lev_files(iff)) THEN 1921 zx_tmp_fi3d(1:klon,1:klev)=d_t_diss(1:klon,1:klev)/pdtphys 1922 CALL histwrite_phy(nid_files(iff),clef_stations(iff), 1923 $o_dtdis%name,itau_w,zx_tmp_fi3d) 1924 ENDIF 1925 1910 1926 IF (o_dqvdf%flag(iff)<=lev_files(iff)) THEN 1911 1927 zx_tmp_fi3d(1:klon,1:klev)=d_q_vdf(1:klon,1:klev)/pdtphys … … 2191 2207 IF (o_evu%flag(iff)<=lev_files(iff)) THEN 2192 2208 CALL histwrite_phy(nid_files(iff),clef_stations(iff), 2193 $o_evu%name,itau_w,coefm )2209 $o_evu%name,itau_w,coefm(:,:,is_ave)) 2194 2210 ENDIF 2195 2211 -
LMDZ5/trunk/libf/phylmd/phys_state_var_mod.F90
r1742 r1761 63 63 !$OMP THREADPRIVATE(ratqs) 64 64 REAL, ALLOCATABLE, SAVE :: pbl_tke(:,:,:) ! turb kinetic energy 65 !$OMP THREADPRIVATE(pbl_tke) 65 REAL, ALLOCATABLE, SAVE :: coefh(:,:,:) ! Kz enthalpie 66 REAL, ALLOCATABLE, SAVE :: coefm(:,:,:) ! Kz momentum 67 !$OMP THREADPRIVATE(pbl_tke, coefh,coefm) 66 68 REAL, ALLOCATABLE, SAVE :: zmax0(:), f0(:) ! 67 69 !$OMP THREADPRIVATE(zmax0,f0) … … 394 396 ALLOCATE(clwcon(klon,klev),rnebcon(klon,klev)) 395 397 ALLOCATE(ratqs(klon,klev)) 396 ALLOCATE(pbl_tke(klon,klev+1,nbsrf)) 398 ALLOCATE(pbl_tke(klon,klev+1,nbsrf+1)) 399 ALLOCATE(coefh(klon,klev+1,nbsrf+1)) 400 ALLOCATE(coefm(klon,klev+1,nbsrf+1)) 397 401 ALLOCATE(zmax0(klon), f0(klon)) 398 402 ALLOCATE(ema_work1(klon,klev), ema_work2(klon,klev)) … … 530 534 deallocate( u_ancien, v_ancien ) 531 535 deallocate( tr_ancien) !RomP 532 deallocate(ratqs, pbl_tke )536 deallocate(ratqs, pbl_tke,coefh,coefm) 533 537 deallocate(zmax0, f0) 534 538 deallocate(ema_work1, ema_work2) -
LMDZ5/trunk/libf/phylmd/physiq.F
r1753 r1761 655 655 real w0(klon) ! Vitesse des thermiques au LCL 656 656 real w_conv(klon) ! Vitesse verticale de grande \'echelle au LCL 657 real tke0(klon,klev+1) ! TKE au début du pas de temps 657 658 real therm_tke_max0(klon) ! TKE dans les thermiques au LCL 658 659 real env_tke_max0(klon) ! TKE dans l'environnement au LCL … … 694 695 cAA 695 696 cAA Pour phytrac 696 cAA697 REAL coefh(klon,klev) ! coef d'echange pour phytrac, valable pour 2<=k<=klev698 REAL coefm(klon,klev) ! coef d'echange pour U, V699 697 REAL u1(klon) ! vents dans la premiere couche U 700 698 REAL v1(klon) ! vents dans la premiere couche V … … 1288 1286 pbase=0 1289 1287 cIM 180608 1290 c pmflxr=0.1291 c pmflxs=0.1292 1288 1293 1289 itau_con=0 … … 1396 1392 1397 1393 CALL phyetat0 ("startphy.nc",clesphy0,tabcntr0) 1394 IF (klon_glo==1) THEN 1395 coefh=0. ; coefm=0. ; pbl_tke=0. 1396 coefh(:,2,:)=1.e-2 ; coefm(:,2,:)=1.e-2 ; pbl_tke(:,2,:)=1.e-2 1397 PRINT*,'FH WARNING : lignes a supprimer' 1398 ENDIF 1398 1399 cIM begin 1399 1400 print*,'physiq: clwcon rnebcon ratqs',clwcon(1,1),rnebcon(1,1) … … 1770 1771 d1a(:,:)=0. 1771 1772 dam(:,:)=0. 1773 pmflxr=0. 1774 pmflxs=0. 1772 1775 ! RomP <<< 1773 1776 … … 1785 1788 ENDDO 1786 1789 ENDDO 1790 tke0(:,:)=pbl_tke(:,:,is_ave) 1787 1791 IF (nqtot.GE.3) THEN 1788 1792 DO iq = 3, nqtot … … 2073 2077 s albsol1, albsol2, sens, evap, 2074 2078 s zxtsol, zxfluxlat, zt2m, qsat2m, 2075 s d_t_vdf, d_q_vdf, d_u_vdf, d_v_vdf, 2079 s d_t_vdf, d_q_vdf, d_u_vdf, d_v_vdf, d_t_diss, 2076 2080 s coefh, coefm, slab_wfbils, 2077 2081 d qsol, zq2m, s_pblh, s_lcl, … … 2089 2093 !----------------------------------------------------------------------------------------- 2090 2094 ! ajout des tendances de la diffusion turbulente 2091 CALL add_phys_tend(d_u_vdf,d_v_vdf,d_t_vdf,d_q_vdf,dql0,'vdf') 2095 CALL add_phys_tend 2096 s (d_u_vdf,d_v_vdf,d_t_vdf+d_t_diss,d_q_vdf,dql0,'vdf') 2092 2097 !----------------------------------------------------------------------------------------- 2093 2098 … … 3175 3180 $ paprs, 3176 3181 $ pplay, 3177 $ coefh ,3182 $ coefh(:,:,is_ave), 3178 3183 $ pphi, 3179 3184 $ t_seri, … … 3696 3701 I paprs, pplay, pmfu, pmfd, 3697 3702 I pen_u, pde_u, pen_d, pde_d, 3698 I cdragh, coefh , fm_therm, entr_therm,3703 I cdragh, coefh(:,:,is_ave), fm_therm, entr_therm, 3699 3704 I u1, v1, ftsol, pctsrf, 3700 3705 I ustar, u10m, v10m, … … 3723 3728 I t,pmfu, pmfd, pen_u, pde_u, pen_d, pde_d, 3724 3729 I fm_therm,entr_therm, 3725 I cdragh,coefh ,u1,v1,ftsol,pctsrf,3730 I cdragh,coefh(:,:,is_ave),u1,v1,ftsol,pctsrf, 3726 3731 I frac_impa, frac_nucl, 3727 3732 I pphis,airephy,dtime,itap, … … 3758 3763 forall (k=1: llm) exner(:, k) = (pplay(:, k)/paprs(:,1))**RKAPPA 3759 3764 CALL ener_conserv(klon,klev,pdtphys,u,v,t,qx(:,:,ivap), 3760 s u_seri,v_seri,t_seri,q_seri, 3765 s u_seri,v_seri,t_seri,q_seri,pbl_tke(:,:,is_ave)-tke0(:,:), 3761 3766 s zmasse,exner,d_t_ec) 3762 3767 t_seri(:,:)=t_seri(:,:)+d_t_ec(:,:) … … 3963 3968 3964 3969 3970 3965 3971 c============================================================= 3966 3972 ! Separation entre thermiques et non thermiques dans les sorties … … 3991 3997 endif 3992 3998 3993 3994 3999 #include "phys_output_write.h" 3995 4000 … … 4070 4075 ENDIF !if callstats 4071 4076 4072 4073 4077 IF (lafin) THEN 4074 4078 itau_phy = itau_phy + itap
Note: See TracChangeset
for help on using the changeset viewer.