Changeset 2901 for LMDZ5/trunk/libf/phylmd
- Timestamp:
- Jun 2, 2017, 6:12:21 PM (7 years ago)
- Location:
- LMDZ5/trunk/libf/phylmd
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/phylmd/cv3_routines.F90
r2818 r2901 111 111 ok_optim_yield=.False. 112 112 CALL getin_p('ok_optim_yield',ok_optim_yield) 113 ok_homo_tend=.TRUE. 114 CALL getin_p('ok_homo_tend',ok_homo_tend) 115 113 116 coef_peel=0.25 114 117 CALL getin_p('coef_peel',coef_peel) … … 2984 2987 qcondc, wd, & 2985 2988 ftd, fqd, qnk, qtc, sigt, tau_cld_cv, coefw_cld_cv) 2989 2990 USE print_control_mod, ONLY: lunout, prt_level 2986 2991 2987 2992 IMPLICIT NONE … … 3873 3878 3874 3879 !!!! do 700 i=1,icb(il)-1 3875 DO i = 1, nl 3876 DO il = 1, ncum 3877 IF (i<=(icb(il)-1) .AND. iflag(il)<=1) THEN 3878 ftd(il, i) = esum(il)*t_wake(il, i)/(th_wake(il,i)*hsum(il)) 3879 fqd(il, i) = fsum(il)/gsum(il) 3880 ft(il, i) = ftd(il, i) + asum(il)*t(il, i)/(th(il,i)*dsum(il)) 3881 fr(il, i) = fqd(il, i) + bsum(il)/csum(il) 3882 END IF 3883 END DO 3884 END DO 3880 IF (ok_homo_tend) THEN 3881 DO i = 1, nl 3882 DO il = 1, ncum 3883 IF (i<=(icb(il)-1) .AND. iflag(il)<=1) THEN 3884 ftd(il, i) = esum(il)*t_wake(il, i)/(th_wake(il,i)*hsum(il)) 3885 fqd(il, i) = fsum(il)/gsum(il) 3886 ft(il, i) = ftd(il, i) + asum(il)*t(il, i)/(th(il,i)*dsum(il)) 3887 fr(il, i) = fqd(il, i) + bsum(il)/csum(il) 3888 END IF 3889 END DO 3890 END DO 3891 ENDIF 3885 3892 3886 3893 !jyg< … … 3920 3927 END DO 3921 3928 ! 3922 ! print *,' YIELD : alpha_qpos ',alpha_qpos(1) 3929 IF (prt_level .GE. 5) THEN 3930 print *,' CV3_YIELD : alpha_qpos ',alpha_qpos(1) 3931 ENDIF 3923 3932 ! 3924 3933 DO il = 1, ncum -
LMDZ5/trunk/libf/phylmd/cv3param.h
r2761 r2901 7 7 !------------------------------------------------------------ 8 8 9 logical ok_homo_tend 9 10 logical ok_optim_yield 10 11 logical ok_convstop … … 41 42 ,noff, minorig, nl, nlp, nlm & 42 43 ,ok_convstop, ok_intermittent & 43 ,ok_optim_yield 44 ,ok_optim_yield & 45 ,ok_homo_tend 44 46 !$OMP THREADPRIVATE(/cv3param/) 45 47
Note: See TracChangeset
for help on using the changeset viewer.