Changeset 3856
- Timestamp:
- Feb 25, 2021, 7:02:16 PM (4 years ago)
- Location:
- LMDZ6/trunk/libf/phylmd
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ6/trunk/libf/phylmd/calcratqs.F90
r2534 r3856 2 2 iflag_ratqs,iflag_con,iflag_cld_th,pdtphys, & 3 3 ratqsbas,ratqshaut,ratqsp0,ratqsdp, & 4 tau_ratqs,fact_cldcon, &4 tau_ratqs,fact_cldcon,wake_s, wake_deltaq, & 5 5 ptconv,ptconvth,clwcon0th, rnebcon0th, & 6 6 paprs,pplay,q_seri,zqsat,fm_therm, & 7 ratqs,ratqsc )7 ratqs,ratqsc,ratqs_inter) 8 8 9 9 implicit none … … 27 27 logical, dimension(klon,klev),intent(in) :: ptconv 28 28 real, dimension(klon,klev),intent(in) :: rnebcon0th,clwcon0th 29 29 real, dimension(klon,klev),intent(in) :: wake_deltaq,wake_s 30 30 ! Output 31 real, dimension(klon,klev),intent(inout) :: ratqs,ratqsc 31 real, dimension(klon,klev),intent(inout) :: ratqs,ratqsc,ratqs_inter 32 32 logical, dimension(klon,klev),intent(inout) :: ptconvth 33 33 … … 124 124 enddo 125 125 126 else if (iflag_ratqs==4) then 126 else if (iflag_ratqs==4) then 127 127 do k=1,klev 128 128 ratqss(:,k)=ratqsbas+0.5*(ratqshaut-ratqsbas) & … … 131 131 enddo 132 132 133 else if (iflag_ratqs==10) then ! ratqs interactif dépendant de la présence de poches froides 134 call calcratqs_inter(klon,klev,pdtphys,ratqsbas,wake_deltaq,wake_s,q_seri,ratqs_inter) 135 do k=1,klev 136 do i=1, klon 137 ratqss(i,k)=ratqs_inter(i,k)+0.5*(ratqshaut-ratqs_inter(i,k)) & 138 *( tanh( (ratqsp0-pplay(i,k))/ratqsdp) + 1.) 139 enddo 140 enddo 141 133 142 endif 134 135 136 143 137 144 -
LMDZ6/trunk/libf/phylmd/phyetat0.F90
r3815 r3856 19 19 wake_s, wake_dens, zgam, zmax0, zmea, zpic, zsig, & 20 20 zstd, zthe, zval, ale_bl, ale_bl_trig, alp_bl, u10m, v10m, treedrg, & 21 ale_wake, ale_bl_stat, ds_ns, dt_ns, delta_sst, delta_sal 21 ale_wake, ale_bl_stat, ds_ns, dt_ns, delta_sst, delta_sal, ratqs_inter 22 22 !FC 23 23 USE geometry_mod, ONLY : longitude_deg, latitude_deg … … 435 435 found=phyetat0_get(1,ale_bl_stat,"ALE_BL_STAT","ALE_BL_STAT",0.) 436 436 437 ! fisrtilp/Clouds 438 found=phyetat0_get(klev,ratqs_inter,"RATQS_INTER","Relative width of the lsc sugrid scale water",ratqs_inter) 439 437 440 !=========================================== 438 441 ! Read and send field trs to traclmdz -
LMDZ6/trunk/libf/phylmd/phyredem.F90
r3815 r3856 28 28 du_gwd_rando, du_gwd_front, u10m, v10m, & 29 29 treedrg, solswfdiff, delta_sal, ds_ns, dt_ns, & 30 delta_sst 30 delta_sst, ratqs_inter 31 31 32 32 USE geometry_mod, ONLY : longitude_deg, latitude_deg … … 303 303 304 304 CALL put_field(pass,"ALE_BL_STAT", "ALE_BL_STAT", ale_bl_stat) 305 306 307 ! fisrtilp/clouds 308 CALL put_field(pass,"RATQS_INTER","Relative width of the lsc sugrid scale water",ratqs_inter) 305 309 306 310 -
LMDZ6/trunk/libf/phylmd/phys_state_var_mod.F90
r3815 r3856 418 418 !$OMP THREADPRIVATE(ccm) 419 419 420 !!! nrlmd le 10/04/2012421 420 REAL,SAVE,ALLOCATABLE :: ale_bl_trig(:) 422 421 !$OMP THREADPRIVATE(ale_bl_trig) 423 !!! fin nrlmd le 10/04/2012 422 423 REAL,SAVE,ALLOCATABLE :: ratqs_inter(:,:) 424 !$OMP THREADPRIVATE(ratqs_inter) 424 425 425 426 REAL, ALLOCATABLE, SAVE:: du_gwd_rando(:, :), du_gwd_front(:, :) … … 648 649 ALLOCATE(cg_aero_lw_rrtm(klon,klev,2,nbands_lw_rrtm)) 649 650 ALLOCATE(ccm(klon,klev,nbands)) 650 651 !!! nrlmd le 10/04/2012652 651 ALLOCATE(ale_bl_trig(klon)) 653 !!! fin nrlmd le 10/04/2012 652 ALLOCATE(ratqs_inter(klon,klev)) 654 653 IF (ok_gwd_rando) THEN 655 654 ALLOCATE(du_gwd_rando(klon, klev)) … … 794 793 if (ok_gwd_rando) DEALLOCATE(du_gwd_rando) 795 794 if (.not. ok_hines .and. ok_gwd_rando) DEALLOCATE(du_gwd_front) 796 797 !!! nrlmd le 10/04/2012798 795 DEALLOCATE(ale_bl_trig) 799 !!! fin nrlmd le 10/04/2012 796 DEALLOCATE(ratqs_inter) 800 797 801 798 if (activate_ocean_skin >= 1) deALLOCATE(delta_sal, ds_ns, dt_ns, & -
LMDZ6/trunk/libf/phylmd/physiq_mod.F90
r3817 r3856 3447 3447 iflag_ratqs,iflag_con,iflag_cld_th,pdtphys, & 3448 3448 ratqsbas,ratqshaut,ratqsp0, ratqsdp, & 3449 tau_ratqs,fact_cldcon, &3449 tau_ratqs,fact_cldcon,wake_s, wake_deltaq, & 3450 3450 ptconv,ptconvth,clwcon0th, rnebcon0th, & 3451 3451 paprs,pplay,q_seri,zqsat,fm_therm, & 3452 ratqs,ratqsc) 3453 3452 ratqs,ratqsc,ratqs_inter) 3454 3453 3455 3454 !
Note: See TracChangeset
for help on using the changeset viewer.