Changeset 2799 for LMDZ5/trunk/libf
- Timestamp:
- Feb 24, 2017, 7:50:33 PM (8 years ago)
- Location:
- LMDZ5/trunk/libf/phylmd
- Files:
-
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/phylmd/YOETHF.h
r2635 r2799 19 19 REAL R5ALVCP,R5ALSCP,RALVDCP,RALSDCP,RALFDCP,RTWAT,RTBER,RTBERCU 20 20 REAL RTICE,RTICECU,RTWAT_RTICE_R,RTWAT_RTICECU_R,RKOOP1,RKOOP2 21 LOGICAL OK_BAD_ECMWF_THERMO ! If TRUE, then variables set by rrtm/suphec.F90 22 ! If FALSE, then variables set by suphel.F90 21 23 COMMON /YOETHF/R2ES, R3LES, R3IES, R4LES, R4IES, R5LES, R5IES, & 22 24 & RVTMP2, RHOH2O, & … … 24 26 & RALFDCP,RTWAT,RTBER,RTBERCU, & 25 27 & RTICE,RTICECU,RTWAT_RTICE_R,RTWAT_RTICECU_R,RKOOP1,& 26 & RKOOP2 28 & RKOOP2, & 29 & OK_BAD_ECMWF_THERMO 27 30 28 31 !$OMP THREADPRIVATE(/YOETHF/) -
LMDZ5/trunk/libf/phylmd/climb_hq_mod.F90
r2159 r2799 9 9 SAVE 10 10 PRIVATE 11 PUBLIC :: climb_hq_down, climb_hq_up 11 PUBLIC :: climb_hq_down, climb_hq_up, d_h_col_vdf, f_h_bnd 12 12 13 13 REAL, DIMENSION(:,:), ALLOCATABLE :: gamaq, gamah … … 23 23 REAL, DIMENSION(:,:), ALLOCATABLE :: Kcoefhq 24 24 !$OMP THREADPRIVATE(Kcoefhq) 25 REAL, SAVE, DIMENSION(:,:), ALLOCATABLE :: h_old ! for diagnostics, h before solving diffusion 26 !$OMP THREADPRIVATE(h_old) 27 REAL, SAVE, DIMENSION(:), ALLOCATABLE :: d_h_col_vdf ! for diagnostics, vertical integral of enthalpy change 28 !$OMP THREADPRIVATE(d_h_col_vdf) 29 REAL, SAVE, DIMENSION(:), ALLOCATABLE :: f_h_bnd ! for diagnostics, enthalpy flux at surface 30 !$OMP THREADPRIVATE(d_h_bnd) 25 31 26 32 CONTAINS … … 71 77 LOGICAL, SAVE :: first=.TRUE. 72 78 !$OMP THREADPRIVATE(first) 73 REAL, DIMENSION(klon,klev) :: local_H 79 ! JLD now renamed h_old and declared in module 80 ! REAL, DIMENSION(klon,klev) :: local_H 74 81 REAL, DIMENSION(klon) :: psref 75 82 REAL :: delz, pkh 76 83 INTEGER :: k, i, ierr 77 78 84 ! Include 79 85 !**************************************************************************************** … … 113 119 ALLOCATE(gamah(1:klon,2:klev), STAT=ierr) 114 120 IF ( ierr /= 0 ) PRINT*,' pb in allloc gamah, ierr=', ierr 121 122 ALLOCATE(h_old(klon,klev), STAT=ierr) 123 IF ( ierr /= 0 ) PRINT*,' pb in allloc h_old, ierr=', ierr 124 125 ALLOCATE(d_h_col_vdf(klon), STAT=ierr) 126 IF ( ierr /= 0 ) PRINT*,' pb in allloc d_h_col_vdf, ierr=', ierr 127 128 ALLOCATE(f_h_bnd(klon), STAT=ierr) 129 IF ( ierr /= 0 ) PRINT*,' pb in allloc f_h_bnd, ierr=', ierr 115 130 END IF 116 131 … … 177 192 ! 178 193 !**************************************************************************************** 179 local_H(:,:) = 0.0194 h_old(:,:) = 0.0 180 195 181 196 DO k=1,klev 182 197 DO i = 1, knon 183 198 ! convertie la temperature en entalpie potentielle 184 local_H(i,k) = RCPD * temp(i,k) * &199 h_old(i,k) = RCPD * temp(i,k) * & 185 200 (psref(i)/pplay(i,k))**RKAPPA 186 201 ENDDO 187 202 ENDDO 188 203 189 CALL calc_coef(knon, Kcoefhq(:,:), gamah(:,:), delp(:,:), local_H(:,:), &204 CALL calc_coef(knon, Kcoefhq(:,:), gamah(:,:), delp(:,:), h_old(:,:), & 190 205 Ccoef_H(:,:), Dcoef_H(:,:), Acoef_H, Bcoef_H) 191 206 … … 349 364 ! 1) 350 365 ! Definition of some variables 366 REAL, DIMENSION(klon,klev) :: d_h, zairm 351 367 ! 352 368 !**************************************************************************************** … … 355 371 d_q(:,:) = 0.0 356 372 d_t(:,:) = 0.0 373 d_h(:,:) = 0.0 374 f_h_bnd(:)= 0.0 357 375 358 376 psref(1:knon) = paprs(1:knon,1) … … 393 411 q_new(1:knon,1) = Acoef_Q(1:knon) + Bcoef_Q(1:knon)*flx_q1(1:knon)*dtime 394 412 h_new(1:knon,1) = Acoef_H(1:knon) + Bcoef_H(1:knon)*flx_h1(1:knon)*dtime 395 413 f_h_bnd(1:knon) = flx_h1(1:knon) 396 414 !- All the other layers 397 415 DO k = 2, klev … … 427 445 ! 428 446 !**************************************************************************************** 429 447 d_h_col_vdf(:) = 0.0 430 448 DO k = 1, klev 431 449 DO i = 1, knon 432 450 d_t(i,k) = h_new(i,k)/(psref(i)/pplay(i,k))**RKAPPA/RCPD - t_old(i,k) 433 451 d_q(i,k) = q_new(i,k) - q_old(i,k) 434 END DO 452 d_h(i,k) = h_new(i,k) - h_old(i,k) 453 !JLD d_t(i,k) = d_h(i,k)/(psref(i)/pplay(i,k))**RKAPPA/RCPD !correction a venir 454 ! layer air mass 455 zairm(i, k) = (paprs(i,k)-paprs(i,k+1))/rg 456 d_h_col_vdf(i) = d_h_col_vdf(i) + d_h(i,k)*zairm(i,k) 457 END DO 435 458 END DO 436 459 … … 448 471 DEALLOCATE(Kcoefhq,stat=ierr) 449 472 IF ( ierr /= 0 ) PRINT*,' pb in dealllocate Kcoefhq, ierr=', ierr 473 DEALLOCATE(h_old, d_h_col_vdf, f_h_bnd, stat=ierr) 474 IF ( ierr /= 0 ) PRINT*,' pb in dealllocate h_old, d_h_col_vdf, f_h_bnd, ierr=', ierr 450 475 END IF 451 476 END SUBROUTINE climb_hq_up -
LMDZ5/trunk/libf/phylmd/phys_output_var_mod.F90
r2752 r2799 24 24 REAL, SAVE, ALLOCATABLE :: bils_latent(:) ! bilan de chaleur au sol 25 25 !$OMP THREADPRIVATE(bils_ec,bils_ech,bils_tke,bils_diss,bils_kinetic,bils_enthalp,bils_latent) 26 ! output variables for energy conservation tests, computed in add_phys_tend 27 REAL, SAVE, ALLOCATABLE :: d_qw_col(:) ! watter vapour mass budget for each column (kg/m2/s) 28 REAL, SAVE, ALLOCATABLE :: d_ql_col(:) ! liquid watter mass budget for each column (kg/m2/s) 29 REAL, SAVE, ALLOCATABLE :: d_qs_col(:) ! solid watter mass budget for each column (kg/m2/s) 30 REAL, SAVE, ALLOCATABLE :: d_qt_col(:) ! total watter mass budget for each column (kg/m2/s) 31 REAL, SAVE, ALLOCATABLE :: d_ek_col(:) ! kinetic energy budget for each column (W/m2) 32 REAL, SAVE, ALLOCATABLE :: d_h_dair_col(:) ! enthalpy budget of dry air for each column (W/m2) 33 REAL, SAVE, ALLOCATABLE :: d_h_qw_col(:) ! enthalpy budget of watter vapour for each column (W/m2) 34 REAL, SAVE, ALLOCATABLE :: d_h_ql_col(:) ! enthalpy budget of liquid watter for each column (W/m2) 35 REAL, SAVE, ALLOCATABLE :: d_h_qs_col(:) ! enthalpy budget of solid watter for each column (W/m2) 36 REAL, SAVE, ALLOCATABLE :: d_h_col(:) ! total enthalpy budget for each column (W/m2) 37 !$OMP THREADPRIVATE(d_qw_col, d_ql_col, d_qs_col, d_qt_col, d_ek_col, d_h_dair_col) 38 !$OMP THREADPRIVATE(d_h_qw_col, d_h_ql_col, d_h_qs_col, d_h_col 26 39 27 40 ! Marine … … 121 134 122 135 allocate (bils_ec(klon),bils_ech(klon),bils_tke(klon),bils_diss(klon),bils_kinetic(klon),bils_enthalp(klon),bils_latent(klon)) 136 allocate (d_qw_col(klon), d_ql_col(klon), d_qs_col(klon), d_qt_col(klon), d_ek_col(klon), d_h_dair_col(klon) & 137 & , d_h_qw_col(klon), d_h_ql_col(klon), d_h_qs_col(klon), d_h_col(klon)) 138 d_qw_col=0. ; d_ql_col=0. ; d_qs_col=0. ; d_qt_col=0. ; d_ek_col=0. ; d_h_dair_col =0. 139 d_h_qw_col=0. ; d_h_ql_col=0. ; d_h_qs_col=0. ; d_h_col=0. 123 140 124 141 ! Marine … … 155 172 deallocate(snow_o,zfra_o,itau_con) 156 173 deallocate (bils_ec,bils_ech,bils_tke,bils_diss,bils_kinetic,bils_enthalp,bils_latent) 174 deallocate (d_qw_col, d_ql_col, d_qs_col, d_qt_col, d_ek_col, d_h_dair_col & 175 & , d_h_qw_col, d_h_ql_col, d_h_qs_col, d_h_col) 157 176 158 177 ! Marine -
LMDZ5/trunk/libf/phylmd/physiq_mod.F90
r2795 r2799 242 242 243 243 USE cmp_seri_mod 244 USE add_phys_tend_mod, only : add_pbl_tend, add_phys_tend, prt_enerbil, & 245 & fl_ebil, fl_cor_ebil 244 246 245 247 !IM stations CFMIP … … 406 408 REAL pphis(klon) 407 409 REAL presnivs(klev) 408 REAL znivsig(klev)409 real pir410 !JLD REAL znivsig(klev) 411 !JLD real pir 410 412 411 413 REAL u(klon,klev) … … 467 469 ! jmin_debut : indice minimum de j; nbptj : nombre de points en 468 470 ! direction j (latitude) 469 INTEGER imin_debut, nbpti470 INTEGER jmin_debut, nbptj471 !JLD INTEGER imin_debut, nbpti 472 !JLD INTEGER jmin_debut, nbptj 471 473 !IM: region='3d' <==> sorties en global 472 474 CHARACTER*3 region … … 617 619 real env_tke_max0(klon) ! TKE dans l'environnement au LCL 618 620 619 !---D\'eclenchement stochastique620 integer :: tau_trig(klon)621 !JLD !---D\'eclenchement stochastique 622 !JLD integer :: tau_trig(klon) 621 623 622 624 REAL,SAVE :: random_notrig_max=1. … … 652 654 REAL beta_prec_fisrt(klon,klev) ! taux de conv de l'eau cond (fisrt) 653 655 ! RomP <<< 654 655 REAL :: calday656 656 657 657 !IM cf FH pour Tiedtke 080604 … … 751 751 REAL conv_t(klon,klev) ! convergence de la temperature(K/s) 752 752 ! 753 #ifdef INCA 753 754 REAL zxsnow_dummy(klon) 755 #endif 754 756 REAL zsav_tsol(klon) 755 757 ! … … 762 764 LOGICAL zx_ajustq 763 765 ! 764 REAL za , zb765 REAL zx_t, zx_qs, zdelta, zcor , zlvdcp, zlsdcp766 REAL za 767 REAL zx_t, zx_qs, zdelta, zcor 766 768 real zqsat(klon,klev) 767 769 ! 768 INTEGER i, k, iq, j, nsrf, ll, l770 INTEGER i, k, iq, nsrf, l 769 771 ! 770 772 REAL t_coup … … 874 876 875 877 ! 876 integer itau_w ! pas de temps ecriture = itap + itau_phy878 !JLD integer itau_w ! pas de temps ecriture = itap + itau_phy 877 879 ! 878 880 ! … … 886 888 REAL tabcntr0( length ) 887 889 ! 888 INTEGER ndex2d(nbp_lon*nbp_lat)890 !JLD INTEGER ndex2d(nbp_lon*nbp_lat) 889 891 !IM 890 892 ! 891 893 !IM AMIP2 BEG 892 REAL moyglo, mountor894 !JLD REAL moyglo, mountor 893 895 !IM 141004 BEG 894 896 REAL zustrdr(klon), zvstrdr(klon) … … 904 906 ! REAL airedyn(nbp_lon+1,nbp_lat) 905 907 !IM 190504 END 906 LOGICAL ok_msk907 REAL msk(klon)908 !JLD LOGICAL ok_msk 909 !JLD REAL msk(klon) 908 910 !ym A voir plus tard 909 911 !ym REAL zm_wo(jjmp1, klev) … … 912 914 REAL zx_tmp_fi2d(klon) ! variable temporaire grille physique 913 915 REAL zx_tmp_fi3d(klon,klev) ! variable temporaire pour champs 3D 914 REAL zx_tmp_2d(nbp_lon,nbp_lat)915 REAL zx_lon(nbp_lon,nbp_lat)916 REAL zx_lat(nbp_lon,nbp_lat)916 !JLD REAL zx_tmp_2d(nbp_lon,nbp_lat) 917 !JLD REAL zx_lon(nbp_lon,nbp_lat) 918 !JLD REAL zx_lat(nbp_lon,nbp_lat) 917 919 ! 918 920 INTEGER nid_ctesGCM … … 930 932 REAL uq_lay(klon,klev) ! transport zonal de l'eau a chaque niveau vert. 931 933 ! 932 REAL zjulian933 SAVE zjulian934 ! $OMP THREADPRIVATE(zjulian)935 936 INTEGER nhori, nvert937 REAL zsto938 REAL zstophy, zout934 !JLD REAL zjulian 935 !JLD SAVE zjulian 936 !JLD!$OMP THREADPRIVATE(zjulian) 937 938 !JLD INTEGER nhori, nvert 939 !JLD REAL zsto 940 !JLD REAL zstophy, zout 939 941 940 942 character*20 modname … … 950 952 parameter (prof2d_on = 1, prof3d_on = 2, & 951 953 prof2d_av = 3, prof3d_av = 4) 952 ! Variables liees au bilan d'energie et d'enthalpi953 954 REAL ztsol(klon) 954 REAL d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec955 REAL d_h_vcol_phy956 REAL fs_bound, fq_bound957 SAVE d_h_vcol_phy958 !$OMP THREADPRIVATE(d_h_vcol_phy)959 REAL zero_v(klon)960 CHARACTER*40 ztit961 INTEGER ip_ebil ! PRINT level for energy conserv. diag.962 SAVE ip_ebil963 DATA ip_ebil/0/964 !$OMP THREADPRIVATE(ip_ebil)965 955 REAL q2m(klon,nbsrf) ! humidite a 2m 966 956 967 957 !IM: t2m, q2m, ustar, u10m, v10m et t2mincels, t2maxcels 968 958 CHARACTER*40 t2mincels, t2maxcels !t2m min., t2m max 969 CHARACTER*40 tinst, tave , typeval959 CHARACTER*40 tinst, tave 970 960 REAL cldtaupi(klon,klev) ! Cloud optical thickness for 971 961 ! pre-industrial (pi) aerosols … … 1025 1015 ! vars_climoz(1:read_climoz): variables names in climoz file. 1026 1016 1027 INTEGER ierr1028 1017 include "YOMCST.h" 1029 1018 include "YOETHF.h" … … 1039 1028 ! Declarations pour Simulateur COSP 1040 1029 !============================================================ 1030 #ifdef CPP_COSP 1041 1031 real :: mr_ozone(klon,klev) 1042 1032 #endif 1043 1033 !IM stations CFMIP 1044 1034 INTEGER, SAVE :: nCFMIP … … 1089 1079 !--OB variables for mass fixer (hard coded for now) 1090 1080 logical, parameter :: mass_fixer=.false. 1091 real qql1(klon),qql2(klon), zdz,corrqql1081 real qql1(klon),qql2(klon),corrqql 1092 1082 1093 1083 ! Ehouarn: set value of jjmp1 since it is no longer a "fixed parameter" … … 1193 1183 1194 1184 modname = 'physiq' 1195 !IM1196 IF (ip_ebil_phy.ge.1) THEN1197 DO i=1,klon1198 zero_v(i)=0.1199 ENDDO1200 ENDIF1201 1185 1202 1186 IF (debut) THEN … … 1210 1194 iflag_wake_tend = 0 1211 1195 CALL getin_p('iflag_wake_tend',iflag_wake_tend) 1196 ok_bad_ecmwf_thermo=.TRUE. ! By default thermodynamical constants are set 1197 ! in rrtm/suphec.F90 (and rvtmp2 is set to 0). 1198 CALL getin_p('ok_bad_ecmwf_thermo',ok_bad_ecmwf_thermo) 1199 fl_ebil = 0 ! by default, conservation diagnostics are desactivated 1200 CALL getin_p('fl_ebil',fl_ebil) 1201 fl_cor_ebil = 0 ! by default, no correction to ensure energy conservation 1202 CALL getin_p('fl_cor_ebil',fl_cor_ebil) 1212 1203 ENDIF 1213 1204 … … 1273 1264 clwcon(:,:) = 0.0 1274 1265 1275 !IM1276 IF (ip_ebil_phy.ge.1) d_h_vcol_phy=0.1277 1266 ! 1278 1267 print*,'iflag_coupl,iflag_clos,iflag_wake', & … … 2008 1997 CALL add_phys_tend & 2009 1998 (du0,dv0,d_t_eva,d_q_eva,d_ql_eva,d_qi_eva,paprs,& 2010 'eva',abortphy,flag_inhib_tend) 1999 'eva',abortphy,flag_inhib_tend,itap) 2000 call prt_enerbil('eva',itap) 2011 2001 2012 2002 !========================================================================= … … 2230 2220 CALL add_pbl_tend & 2231 2221 (d_u_vdf,d_v_vdf,d_t_vdf+d_t_diss,d_q_vdf,dql0,dqi0,paprs,& 2232 'vdf',abortphy,flag_inhib_tend )2222 'vdf',abortphy,flag_inhib_tend,itap) 2233 2223 ELSE 2234 2224 CALL add_phys_tend & 2235 2225 (d_u_vdf,d_v_vdf,d_t_vdf+d_t_diss,d_q_vdf,dql0,dqi0,paprs,& 2236 'vdf',abortphy,flag_inhib_tend) 2237 ENDIF 2226 'vdf',abortphy,flag_inhib_tend,itap) 2227 ENDIF 2228 call prt_enerbil('vdf',itap) 2238 2229 !-------------------------------------------------------------------- 2239 2230 … … 2616 2607 2617 2608 CALL add_phys_tend(d_u_con, d_v_con, d_t_con, d_q_con, dql0, dqi0, paprs, & 2618 'convection',abortphy,flag_inhib_tend) 2609 'convection',abortphy,flag_inhib_tend,itap) 2610 call prt_enerbil('convection',itap) 2619 2611 2620 2612 !------------------------------------------------------------------------- … … 2745 2737 ! ajout des tendances des poches froides 2746 2738 CALL add_phys_tend(du0,dv0,d_t_wake,d_q_wake,dql0,dqi0,paprs,'wake', & 2747 abortphy,flag_inhib_tend) 2739 abortphy,flag_inhib_tend,itap) 2740 call prt_enerbil('wake',itap) 2748 2741 !------------------------------------------------------------------------ 2749 2742 … … 2754 2747 (d_deltat_wk, d_deltaq_wk, d_s_wk, d_dens_wk, wake_k, & 2755 2748 'wake', abortphy) 2756 2749 call prt_enerbil('wake',itap) 2757 2750 ENDIF ! (iflag_wake_tend .GT. 0.) 2758 2751 … … 2875 2868 CALL add_wake_tend & 2876 2869 (d_deltat_the, d_deltaq_the, dsig0, ddens0, wkoccur1, 'the', abortphy) 2870 call prt_enerbil('the',itap) 2877 2871 ! 2878 2872 ENDIF ! (mod(iflag_pbl_split/2,2) .EQ. 1) 2879 2873 ! 2880 2874 CALL add_phys_tend(d_u_ajs,d_v_ajs,d_t_ajs,d_q_ajs, & 2881 dql0,dqi0,paprs,'thermals', abortphy,flag_inhib_tend) 2875 dql0,dqi0,paprs,'thermals', abortphy,flag_inhib_tend,itap) 2876 call prt_enerbil('thermals',itap) 2882 2877 ! 2883 2878 ! … … 2941 2936 ! ajout des tendances de l'ajustement sec ou des thermiques 2942 2937 CALL add_phys_tend(du0,dv0,d_t_ajsb,d_q_ajsb,dql0,dqi0,paprs, & 2943 'ajsb',abortphy,flag_inhib_tend) 2938 'ajsb',abortphy,flag_inhib_tend,itap) 2939 call prt_enerbil('ajsb',itap) 2944 2940 d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_ajsb(:,:) 2945 2941 d_q_ajs(:,:)=d_q_ajs(:,:)+d_q_ajsb(:,:) … … 2984 2980 WHERE (snow_lsc < 0) snow_lsc = 0. 2985 2981 2982 !+JLD 2983 ! write(*,9000) 'phys lsc',"enerbil: bil_q, bil_e,",rain_lsc+snow_lsc & 2984 ! & ,((rcw-rcpd)*rain_lsc + (rcs-rcpd)*snow_lsc)*t_seri(1,1)-rlvtt*rain_lsc+rlstt*snow_lsc & 2985 ! & ,rain_lsc,snow_lsc 2986 ! write(*,9000) "rcpv","rcw",rcpv,rcw,rcs,t_seri(1,1) 2987 !-JLD 2986 2988 CALL add_phys_tend(du0,dv0,d_t_lsc,d_q_lsc,d_ql_lsc,d_qi_lsc,paprs, & 2987 'lsc',abortphy,flag_inhib_tend) 2989 'lsc',abortphy,flag_inhib_tend,itap) 2990 call prt_enerbil('lsc',itap) 2988 2991 rain_num(:)=0. 2989 2992 DO k = 1, klev … … 3764 3767 ENDDO 3765 3768 3766 CALL add_phys_tend(du0,dv0,d_t_swr,dq0,dql0,dqi0,paprs,'SW',abortphy,flag_inhib_tend) 3767 CALL add_phys_tend(du0,dv0,d_t_lwr,dq0,dql0,dqi0,paprs,'LW',abortphy,flag_inhib_tend) 3769 CALL add_phys_tend(du0,dv0,d_t_swr,dq0,dql0,dqi0,paprs,'SW',abortphy,flag_inhib_tend,itap) 3770 call prt_enerbil('SW',itap) 3771 CALL add_phys_tend(du0,dv0,d_t_lwr,dq0,dql0,dqi0,paprs,'LW',abortphy,flag_inhib_tend,itap) 3772 call prt_enerbil('LW',itap) 3768 3773 3769 3774 ! … … 3835 3840 ! ajout des tendances de la trainee de l'orographie 3836 3841 CALL add_phys_tend(d_u_oro,d_v_oro,d_t_oro,dq0,dql0,dqi0,paprs,'oro', & 3837 abortphy,flag_inhib_tend) 3842 abortphy,flag_inhib_tend,itap) 3843 call prt_enerbil('oro',itap) 3838 3844 !---------------------------------------------------------------------- 3839 3845 ! … … 3881 3887 ! ajout des tendances de la portance de l'orographie 3882 3888 CALL add_phys_tend(d_u_lif, d_v_lif, d_t_lif, dq0, dql0, dqi0, paprs, & 3883 'lif', abortphy,flag_inhib_tend) 3889 'lif', abortphy,flag_inhib_tend,itap) 3890 call prt_enerbil('lif',itap) 3884 3891 ENDIF ! fin de test sur ok_orolf 3885 3892 … … 3904 3911 d_t_hin(:, :)=0. 3905 3912 CALL add_phys_tend(du_gwd_hines, dv_gwd_hines, d_t_hin, dq0, dql0, & 3906 dqi0, paprs, 'hin', abortphy,flag_inhib_tend) 3913 dqi0, paprs, 'hin', abortphy,flag_inhib_tend,itap) 3914 call prt_enerbil('hin',itap) 3907 3915 ENDIF 3908 3916 … … 3921 3929 3922 3930 CALL add_phys_tend(du_gwd_front, dv_gwd_front, dt0, dq0, dql0, dqi0, & 3923 paprs, 'front_gwd_rando', abortphy,flag_inhib_tend) 3931 paprs, 'front_gwd_rando', abortphy,flag_inhib_tend,itap) 3932 call prt_enerbil('front_gwd_rando',itap) 3924 3933 ENDIF 3925 3934 … … 3929 3938 du_gwd_rando, dv_gwd_rando, east_gwstress, west_gwstress) 3930 3939 CALL add_phys_tend(du_gwd_rando, dv_gwd_rando, dt0, dq0, dql0, dqi0, & 3931 paprs, 'flott_gwd_rando', abortphy,flag_inhib_tend) 3940 paprs, 'flott_gwd_rando', abortphy,flag_inhib_tend,itap) 3941 call prt_enerbil('flott_gwd_rando',itap) 3932 3942 zustr_gwd_rando=0. 3933 3943 zvstr_gwd_rando=0. … … 3979 3989 ! ajout de la tendance d'humidite due au methane 3980 3990 CALL add_phys_tend(du0, dv0, dt0, d_q_ch4*dtime, dql0, dqi0, paprs, & 3981 'q_ch4', abortphy,flag_inhib_tend )3991 'q_ch4', abortphy,flag_inhib_tend,itap) 3982 3992 ENDIF 3983 3993 ! -
LMDZ5/trunk/libf/phylmd/reevap.F90
r2705 r2799 2 2 & d_t_eva,d_q_eva,d_ql_eva,d_qs_eva) 3 3 4 4 ! flag to include modifications to ensure energy conservation (if flag >0) 5 USE add_phys_tend_mod, only : fl_cor_ebil 6 5 7 IMPLICIT none 6 8 !>====================================================================== … … 25 27 DO k = 1, klev ! re-evaporation de l'eau liquide nuageuse 26 28 DO i = 1, klon 29 if (fl_cor_ebil .GT. 0) then 30 zlvdcp=RLVTT/RCPD/(1.0+RVTMP2*(q_seri(i,k)+ql_seri(i,k)+qs_seri(i,k))) 31 zlsdcp=RLSTT/RCPD/(1.0+RVTMP2*(q_seri(i,k)+ql_seri(i,k)+qs_seri(i,k))) 32 else 27 33 zlvdcp=RLVTT/RCPD/(1.0+RVTMP2*q_seri(i,k)) 28 34 !jyg< … … 30 36 ! A verifier !!! 31 37 zlsdcp=RLSTT/RCPD/(1.0+RVTMP2*q_seri(i,k)) 38 end if 32 39 IF (iflag_ice_thermo .EQ. 0) THEN 33 40 zlsdcp=zlvdcp -
LMDZ5/trunk/libf/phylmd/rrtm/suphec.F90
r2345 r2799 129 129 130 130 IF (LHOOK) CALL DR_HOOK('SUPHEC',0,ZHOOK_HANDLE) 131 !CALL GSTATS(1811,0) ! MPL 28.11.08 132 !RVTMP2=RCPV/RCPD-1.0_JPRB !use cp,moist 133 RVTMP2=0.0_JPRB !neglect cp,moist 134 RHOH2O=RATM/100._JPRB 135 R2ES=611.21_JPRB*RD/RV 136 R3LES=17.502_JPRB 137 R3IES=22.587_JPRB 138 R4LES=32.19_JPRB 139 R4IES=-0.7_JPRB 140 R5LES=R3LES*(RTT-R4LES) 141 R5IES=R3IES*(RTT-R4IES) 142 R5ALVCP=R5LES*RLVTT/RCPD 143 R5ALSCP=R5IES*RLSTT/RCPD 144 RALVDCP=RLVTT/RCPD 145 RALSDCP=RLSTT/RCPD 146 RALFDCP=RLMLT/RCPD 147 RTWAT=RTT 148 RTBER=RTT-5._JPRB 149 RTBERCU=RTT-5.0_JPRB 150 RTICE=RTT-23._JPRB 151 RTICECU=RTT-23._JPRB 152 153 RTWAT_RTICE_R=1.0_JPRB/(RTWAT-RTICE) 154 RTWAT_RTICECU_R=1.0_JPRB/(RTWAT-RTICECU) 155 IF(NPHYINT == 0) THEN 156 ISMAX=NSMAX 157 ELSE 158 ISMAX=PHYS_GRID%NSMAX 159 ENDIF 160 161 RKOOP1=2.583_JPRB 162 RKOOP2=0.48116E-2_JPRB 131 ! 132 IF (OK_BAD_ECMWF_THERMO) THEN 133 ! 134 ! Modify constants defined in suphel.F90 and set RVTMP2 to 0. 135 ! CALL GSTATS(1811,0) ! MPL 28.11.08 136 ! RVTMP2=RCPV/RCPD-1.0_JPRB !use cp,moist 137 RVTMP2=0.0_JPRB !neglect cp,moist 138 RHOH2O=RATM/100._JPRB 139 R2ES=611.21_JPRB*RD/RV 140 R3LES=17.502_JPRB 141 R3IES=22.587_JPRB 142 R4LES=32.19_JPRB 143 R4IES=-0.7_JPRB 144 R5LES=R3LES*(RTT-R4LES) 145 R5IES=R3IES*(RTT-R4IES) 146 R5ALVCP=R5LES*RLVTT/RCPD 147 R5ALSCP=R5IES*RLSTT/RCPD 148 RALVDCP=RLVTT/RCPD 149 RALSDCP=RLSTT/RCPD 150 RALFDCP=RLMLT/RCPD 151 RTWAT=RTT 152 RTBER=RTT-5._JPRB 153 RTBERCU=RTT-5.0_JPRB 154 RTICE=RTT-23._JPRB 155 RTICECU=RTT-23._JPRB 156 157 RTWAT_RTICE_R=1.0_JPRB/(RTWAT-RTICE) 158 RTWAT_RTICECU_R=1.0_JPRB/(RTWAT-RTICECU) 159 IF(NPHYINT == 0) THEN 160 ISMAX=NSMAX 161 ELSE 162 ISMAX=PHYS_GRID%NSMAX 163 ENDIF 164 165 RKOOP1=2.583_JPRB 166 RKOOP2=0.48116E-2_JPRB 167 168 ELSE 169 ! Keep constants defined in suphel.F90 170 RTICE=RTT-23._JPRB 171 ! 172 ENDIF ! (OK_BAD_ECMWF_THERMO) 163 173 164 174 ! ------------------------------------------------------------------
Note: See TracChangeset
for help on using the changeset viewer.