module process_cat ! reichle, 13 Jun 2005 ! qliu+reichle, 14 Aug 2008 - major overhaul, use subroutine louissurface() ! reichle, 29 Nov 2010 - added Helfand Monin-Obukhov surface layer turbulence scheme ! reichle, 20 Dec 2011 - reinstated "PARdrct" and "PARdffs" for MERRA-Land file specs ! reichle, 28 Dec 2011 - removed field "totalb" from "cat_diagn" structure ! reichle, 15 Mar 2013 - minor changes to adjust to infrastructure changes ! in MAPL_Base.F90 of GEOSadas-5_9_1_p8 ! reichle, 30 Oct 2013 - removed field rzeq from cat_diagn structure ! reichle, 31 Oct 2013 - split "cat_diagn" structure into "cat_diagS" and "cat_diagF" ! reichle, 8 Aug 2014 - use Ganymed-4_0 version of catchment.F90 and StieglitzSnow.F90 ! reichle, 5 Dec 2014 - use Ganymed-4_1 revisions of roughness length use MAPL_SatVaporMod, ONLY: & MAPL_EQsat use catch_types, ONLY: & cat_param_type, & cat_progn_type, & cat_diagF_type, & cat_diagS_type, & catprogn2wesn, & assignment (=) use driver_types, ONLY: & met_force_type, & veg_param_type, & alb_param_type use catchment_model, ONLY: & catchment, & sibalb, & catch_calc_soil_moist, & subtile2tile => catch_calc_subtile2tile use catch_constants, ONLY : & N_snow => CATCH_N_SNOW, & N_gt => CATCH_N_GT, & RHOFS => CATCH_SNWALB_RHOFS, & SNWALB_VISMAX => CATCH_SNWALB_VISMAX, & SNWALB_NIRMAX => CATCH_SNWALB_NIRMAX, & SLOPE => CATCH_SNWALB_SLOPE, & MAPL_NumVegTypes, & Tzero => MAPL_TICE, & stefan_boltzmann => MAPL_STFBOL, & rgas => MAPL_rgas, & grav => MAPL_GRAV, & cpair => MAPL_CP, & VIREPS => MAPL_VIREPS, & MAPL_USMIN, & MAPL_PI, & nodata_generic, & nodata_tol_generic, & sfc_turb_scheme, & GEOS5_coupled_style, & incl_Louis_ww, & incl_Louis_extra_derivs use StieglitzSnow, ONLY: & snowrt, & snow_albedo, & StieglitzSnow_calc_asnow, & StieglitzSnow_calc_tpsnow use sfclayer, ONLY: & louissurface, & helfsurface implicit none private public :: propagate_cat ! parameters for computing zthick, taken from GEOS_CatchGridComp.F90 real, parameter :: MIN_VEG_HEIGHT = 0.01 real, parameter :: SCALE4Z0 = 2. ! default from Ganymed-4_1 real, parameter :: Z0_BY_ZVEG = 0.13 real, parameter :: D0_BY_ZVEG = 0.66 real, parameter :: HPBL = 1000. real, parameter :: DZE_MIN = 10. real, parameter :: TOTALB_MIN = 0.01 real, parameter :: TOTALB_MAX = 0.90 ! parameters for computing land surface emissivity, taken from GEOS_CatchGridComp.F90 ! Emissivity values from Wilber et al (1999, NATA-TP-1999-209362) ! Fu-Liou bands have been combined to Chou bands (though these are broadband only) ! IGBP veg types have been mapped to Sib-Mosaic types ! Details in ~suarez/Emiss on cerebus real, parameter :: EMSVEG(MAPL_NumVegTypes) = (/ & 0.99560, 0.99000, 0.99560, 0.99320, & 0.99280, 0.99180 /) real, parameter :: EMSBARESOIL = 0.94120 real, parameter :: EMSSNO = 0.99999 ! definition of land in sfclayer module integer, parameter :: SFCTYPE = 3 ! parameters for the Helfand Monin-Obuhkhov surface layer routine integer, parameter :: HELF_NITER = 6 ! number of internal iterations integer, parameter :: HELF_Z0 = 1 ! CHOOSEZ0 in GEOS_CatchGridComp.F90 (ocean only) contains subroutine propagate_cat( & N_cat, dtstep, tile_id, & cat_param, & met_force, grn, lai, sunang, alb_param, & ! ntp values cat_progn, & cat_diagS, cat_diagF ) implicit none integer, intent(in) :: N_cat integer, dimension(N_cat), intent(in) :: tile_id real, intent(in) :: dtstep type(cat_param_type), dimension(N_cat), intent(in) :: cat_param type(met_force_type), dimension(N_cat), intent(in) :: met_force !mx type(veg_param_type), dimension(N_cat), intent(in) :: veg_param real, dimension(N_cat), intent(in) :: grn, lai, sunang type(alb_param_type), dimension(N_cat), intent(in) :: alb_param type(cat_progn_type), dimension(N_cat), intent(inout) :: cat_progn type(cat_diagS_type), dimension(N_cat), intent(out) :: cat_diagS type(cat_diagF_type), dimension(N_cat), intent(out) :: cat_diagF ! -------------------------------------------------------------- real :: SFRAC real, dimension(N_cat) :: & SQSCAT, Z2, DZM, RSOIL1, RSOIL2, SATCAP, RDC, U2FAC, & TRAINL, PARDIR, PARDIF, & SWNETF, SWNETS, & RA1, ETURB1, DEDTC1, HSTURB1, DHSDTC1, & RA2, ETURB2, DEDTC2, HSTURB2, DHSDTC2, & RA4, ETURB4, DEDTC4, HSTURB4, DHSDTC4, & RAS, ETURBS, DEDTCS, HSTURBS, DHSDTCS, & QSAT1, DQS1, ALW1, BLW1, & QSAT2, DQS2, ALW2, BLW2, & QSAT4, DQS4, ALW4, BLW4, & QSATS, DQSS, ALWS, BLWS, & DEDQA1, DHSDQA1, & DEDQA2, DHSDQA2, & DEDQA4, DHSDQA4, & DEDQAS, DHSDQAS, & ALBVR, ALBNR, ALBVF, ALBNF, & SNOVR, SNONR, SNOVF, SNONF, & SWLAND, RZEQ, & GHFLXSNO, GHFLXTSKIN, & ENTOT, WTOT, WCHANGE, ECHANGE, & SHSNOW1, AVETSNOW1, WAT10CM1, WATSOI1, ICESOI1, & LHSNOW1, LWUPSNOW1, LWDNSNOW1, NETSWSNOW, & TCSORIG1, TPSN1IN1, TPSN1OUT1, & TC1_0, TC2_0, TC4_0, QA1_0, QA2_0, QA4_0, EACC_0 real, dimension(N_cat) :: zol, zol_sc, d0, dze, tmprealvec, uuu, zvg, ar4, zol_uuu real, dimension(:), allocatable :: ch_tile, cq_tile, qa_tile real, dimension(:,:), allocatable :: qa real, dimension(N_gt, N_cat) :: GHTCNT real, dimension(N_snow,N_cat) :: WESN, HTSN, SNDZ real :: albedo, albsn, asnow, totalb real :: qsattc, tmpt, tmpswnet, emis integer :: k, n logical :: snow_dominates_alb character(len=*), parameter :: Iam = 'propagate_cat' character(len=400) :: err_msg ! ---------------------------------------------------------------- ! initialize cat_diagS and cat_diagF fields that are needed here ! and in sfc_turb !mx ques do n=1,N_cat !mx ques !mx ques cat_diagS(n) = nodata_generic !mx ques cat_diagF(n) = nodata_generic !mx ques !mx ques end do ! diagnose top layer snow temperature in K call StieglitzSnow_calc_tpsnow( N_cat, cat_progn%htsn(1), cat_progn%wesn(1), & cat_diagS%tpsn(1), tmprealvec ) cat_diagS%tpsn(1) = cat_diagS%tpsn(1) + Tzero ! convert to Kelvin ! diagnose sub-tile area fractions ! updated to new interfaces - reichle, 3 Apr 2012 call catch_calc_soil_moist( & N_cat, & cat_param%vegcls, cat_param%dzsf, cat_param%vgwmax, & cat_param%cdcr1, cat_param%cdcr2, cat_param%psis, & cat_param%bee, cat_param%poros, cat_param%wpwet, & cat_param%ars1, cat_param%ars2, cat_param%ars3, & cat_param%ara1, cat_param%ara2, cat_param%ara3, & cat_param%ara4, cat_param%arw1, cat_param%arw2, & cat_param%arw3, cat_param%arw4, & cat_progn%srfexc, cat_progn%rzexc, cat_progn%catdef, & cat_diagS%ar1, cat_diagS%ar2, ar4 ) call StieglitzSnow_calc_asnow( N_snow, N_cat, catprogn2wesn(N_cat,cat_progn), & cat_diagS%asnow ) ! --------------------------------------------------------------- ! assorted Catchment model inputs call pmonth( & N_cat, cat_param%vegcls, grn, lai, & zol, d0, SQSCAT, Z2, zvg, RSOIL1, RSOIL2, RDC ) ! --------------------------------------------------------------- ! precipitation TRAINL = met_force%Rainf - met_force%Rainf_C SFRAC = 1. ! --------------------------------------------------------------- ! turbulence do n=1,N_cat ! Ganymed-4_1 revision of roughness length. Intentionally applies only to ! roughness that is input into surface turbulence calculations, and not to ! wind-related input into subroutine catchment(), see below. ! ! NOTE: Changes in this subroutine (propagate_cat()) to accommodate the ! clean-up of subroutine pmonth() result in roundoff differences ! w.r.t. previous LDASsa version even if SCALE4Z0=1. ! ! - reichle, 8 Dec 2014 zol_sc(n) = zol(n)*SCALE4Z0 zol_uuu(n) = zol(n) dze(n) = max(met_force(n)%RefH-d0(n), DZE_MIN) end do ! compute surface exchange coefficients etc BEFORE possibly resetting ! profile of Qair-QAx-Qsat(surf) -- for consistency with two-stage ! run-method in GEOS_CatchGridComp.F90 ! reichle+qliu, 9 Oct 2008 select case (sfc_turb_scheme) case (0) ! Louis allocate(ch_tile(N_cat)) allocate(cq_tile(N_cat)) call sfc_turb(N_cat, met_force, cat_progn, cat_diagS, zol_sc, dze, lai, & RA1, ETURB1, DEDQA1, DEDTC1, HSTURB1, DHSDQA1, DHSDTC1, & RA2, ETURB2, DEDQA2, DEDTC2, HSTURB2, DHSDQA2, DHSDTC2, & RA4, ETURB4, DEDQA4, DEDTC4, HSTURB4, DHSDQA4, DHSDTC4, & RAS, ETURBS, DEDQAS, DEDTCS, HSTURBS, DHSDQAS, DHSDTCS, & ch_tile=ch_tile, cq_tile=cq_tile ) case (1) ! Helfand Monin-Obukhov call sfc_turb(N_cat, met_force, cat_progn, cat_diagS, zol_sc, dze, lai, & RA1, ETURB1, DEDQA1, DEDTC1, HSTURB1, DHSDQA1, DHSDTC1, & RA2, ETURB2, DEDQA2, DEDTC2, HSTURB2, DHSDQA2, DHSDTC2, & RA4, ETURB4, DEDQA4, DEDTC4, HSTURB4, DHSDQA4, DHSDTC4, & RAS, ETURBS, DEDQAS, DEDTCS, HSTURBS, DHSDQAS, DHSDTCS, & t2m_tile=cat_diagF%t2m, q2m_tile=cat_diagF%q2m ) case default err_msg = 'Unknown selection for sfc_turb_scheme.' stop end select do n=1,N_cat ! make sure that Qair-QAx-Qsat(surf) profile is ok ! (in a nutshell, need Qair<=QAx<=Qsat(surf) or Qair>=QAx>=Qsat(surf)) ! - reichle, 5 Mar 2008 ! - qliu+reichle, 12 Aug 2008 !qsattc = qsat(cat_progn(n)%tc1,PSUR(n),ALHE) qsattc = MAPL_EQsat(cat_progn(n)%tc1,met_force(n)%Psurf) if (cat_progn(n)%qa1>qsattc .and. cat_progn(n)%qa1>met_force(n)%Qair ) then cat_progn(n)%qa1 = max( met_force(n)%Qair, qsattc ) endif if (cat_progn(n)%qa1qsattc .and. cat_progn(n)%qa2>met_force(n)%Qair ) then cat_progn(n)%qa2 = max( met_force(n)%Qair, qsattc ) endif if (cat_progn(n)%qa2qsattc .and. cat_progn(n)%qa4>met_force(n)%Qair ) then cat_progn(n)%qa4 = max( met_force(n)%Qair, qsattc ) endif if (cat_progn(n)%qa40. .and. sunang(n)>0.) then ! sun is up and SWdown>0 ! determine whether snow or snow-free albedo dominates snow_dominates_alb = ( (1-asnow)*(1-albedo)<=asnow*(1-albsn) ) ! compute SWNET for snow-free and snow conditions if ( (abs(met_force(n)%SWnet-nodata_generic) this produces a ! minus sign ! ! taken together, we need d(cq)/dq = d(cq)/dTV * dTV/dQ = -dcq * VIREPS*T dcqdq = -dcq(k,:)*VIREPS*tc(k,:) DEDQA1(k) = cq(k,1) + max( 0., dcqdq(1)*dq(1) ) DEDQA2(k) = cq(k,2) + max( 0., dcqdq(2)*dq(2) ) DEDQA4(k) = cq(k,3) + max( 0., dcqdq(3)*dq(3) ) DEDQAS(k) = cq(k,4) + max( 0., dcqdq(4)*dq(4) ) ! similarly, d(cq)/dt = d(cq)/dTV * dTV/dT = -dcq * (1+VIREPS*Q) dcqdt = -dcq(k,:)*(1.+VIREPS*qa(k,:)) DEDTC1(k) = max( 0., dcqdt(1)*dq(1) ) DEDTC2(k) = max( 0., dcqdt(2)*dq(2) ) DEDTC4(k) = max( 0., dcqdt(3)*dq(3) ) DEDTCS(k) = max( 0., dcqdt(4)*dq(4) ) ! "dch" from louissurface() is d(ch)/d(TVA) where TVA = T*(1+VIREPS*Q) ! ! moreover, ch=ch(Ri) where Ri is proportional to deltaTVA=TVA-TVS ! (where TVS is the virtual surface temperature) -> this produces a ! minus sign ! ! taken together, we need d(ch)/dT = d(ch)/dTV * dTV/dT = -dch * (1+VIREPS*Q) dchdt = -dch(k,:)*(1.+VIREPS*qa(k,:)) DHSDTC1(k) = CPAIR*(ch(k,1) + max( 0., dchdt(1)*dt(1)) ) DHSDTC2(k) = CPAIR*(ch(k,2) + max( 0., dchdt(2)*dt(2)) ) DHSDTC4(k) = CPAIR*(ch(k,3) + max( 0., dchdt(3)*dt(3)) ) DHSDTCS(k) = CPAIR*(ch(k,4) + max( 0., dchdt(4)*dt(4)) ) ! similarly, d(ch)/dQ = d(ch)/dTV * dTV/dQ = -dch * VIREPS*T dchdq = -dch(k,:)*VIREPS*tc(k,:) DHSDQA1(k) = max( 0., CPAIR*dchdq(1)*dt(1) ) DHSDQA2(k) = max( 0., CPAIR*dchdq(2)*dt(2) ) DHSDQA4(k) = max( 0., CPAIR*dchdq(3)*dt(3) ) DHSDQAS(k) = max( 0., CPAIR*dchdq(4)*dt(4) ) end if end do ! T2m, Q2m diagnostics if (sfc_turb_scheme==0 .and. present(ch_tile) .and. present(cq_tile )) then ! Louis: compute tile-avg ch and cq so that T2m, Q2m can be diagnosed later call subtile2tile(N_cat,cat_diagS%ar1,cat_diagS%ar2,cat_diagS%asnow,ch,ch_tile) call subtile2tile(N_cat,cat_diagS%ar1,cat_diagS%ar2,cat_diagS%asnow,cq,cq_tile) end if if (sfc_turb_scheme==1 .and. present(t2m_tile) .and. present(q2m_tile)) then ! Helfand-MO: compute tile-avg T2m, Q2m call subtile2tile(N_cat,cat_diagS%ar1,cat_diagS%ar2,cat_diagS%asnow,t2m,t2m_tile) call subtile2tile(N_cat,cat_diagS%ar1,cat_diagS%ar2,cat_diagS%asnow,q2m,q2m_tile) end if ! ----------------------------------------------------------------- ! ! deallocate select case (sfc_turb_scheme) case (0) ! Louis deallocate(cn) deallocate(ri) deallocate(zt) deallocate(zq) deallocate(uuu) deallocate(ucn) deallocate(re) deallocate(cm) deallocate(ww) deallocate(dch) deallocate(dcq) case (1) ! Helfand Monin-Obukhov deallocate(tmpintvec) deallocate(tmprealvec) deallocate(psmb) deallocate(psl) deallocate(RHOH) deallocate(VKM) deallocate(USTAR) deallocate(XX) deallocate(YY) deallocate(CU) deallocate(CT) deallocate(RIB) deallocate(ZETA) deallocate(WS) deallocate(t2m) deallocate(q2m) deallocate(u2m) deallocate(v2m) deallocate(t10m) deallocate(q10m) deallocate(u10m) deallocate(v10m) deallocate(u50m) deallocate(v50m) end select end subroutine sfc_turb ! ***************************************************************************** SUBROUTINE PMONTH ( & NCH, ITYP, & GREEN, ZLT, Z0, D, & SQSCAT, Z2, ZVG, & RSOIL1, RSOIL2, RDC & ) ! reichle+qliu, 8 Oct 2008 - updated to match GEOS5-MERRA parameters ! reichle 10 Feb 2011 - moved here from file pmonth_stage4.F90 ! reichle 8 Dec 2014 - major clean-up, roundoff differences for two reasons: ! 1) use of MAPL_PI here and ! 2) associated Dec 5 changes in subroutine propagate_cat() IMPLICIT NONE ! This subroutine sets seasonally varying ! variables. NOTE: LAI is stored as LAI (SIB) / VCOVER (SIB). ! INTEGER ntyps parameter (NTYPS = 10) ! ! IntegerTYPe: Vegetation types (biomes). ! 1: BrdEvr: Broadleaf evergreen. ! 2: BrdDcd: Broadleaf deciduous. ! 3: Needle: Needleleaf trees. ! 4: GndCvr: Groundcover (grass). INTEGER, intent(in) :: NCH INTEGER, dimension(NCH), intent(in) :: ITYP REAL, dimension(NCH), intent(in) :: GREEN, ZLT REAL, dimension(NCH), intent(out) :: Z0, D, SQSCAT, Z2, ZVG, RSOIL1, RSOIL2, RDC ! ! local variables ! INTEGER VEG, ChNo ! REAL VGRF11(NTYPS), VGRF12(NTYPS),& VGTR11(NTYPS), VGTR12(NTYPS),& VGROCA(NTYPS), VGROTD(NTYPS),& VGRDRS(NTYPS), VGZ2(NTYPS), & VGRT(NTYPS) ! REAL ALPHAF, DUM1, DUM2, ROOTL, SCAT, VROOT ! ! ----- Sarith 4/19/06 -- from GEOScatchGridComp.F90 ! REAL, DIMENSION (NTYPS-2) :: VGRDA, VGRDB ! Correction to RDC formulation -Randy Koster, 4/1/2011 ! !data VGRDA / 285.9, 294.9, 652.9, 25.8, 100.7, & ! 22.9, 23.8, 23.8/ !data VGRDB / 5.1 , 7.2, 10.8, 4.8, 1.8, 5.1, .000, .000/ ! data VGRDA / 285.9, 355.18, 660.24, 30.06, 100.7, 24.36, 23.8, 23.8/ data VGRDB / 5.1 , 7.2, 10.5, 4.8, 1.8, 5.1, .000, .000/ DATA VGRF11 /0.10, 0.10, 0.07, 0.105 & , 0.10, 0.10, .001, .001, .001, .001/ DATA VGRF12 /0.16, 0.16, 0.16, 0.360& , 0.16, 0.16, .001, .001, .001, .001/ DATA VGTR11 /0.05, 0.05, 0.05, 0.070& , 0.05, 0.05, .001, .001, .001, .001/ DATA VGTR12 /.001, .001, .001, .220& , .001, .001, .001, .001, .001, .001/ DATA VGROCA / & 0.384E-6, 0.384E-6, 0.384E-6, 0.384E-6, 0.384E-6, 0.384E-6,& .1E-6, .1E-6, .1E-6, .1E-6 / DATA VGROTD /1.00, 1.00, 0.50, 0.50, 0.50, 0.20, 0.10, 0.10,& 0.10, 0.10/ DATA VGRDRS /& 0.75E13, 0.75E13, 0.75E13, 0.40E13, 0.75E13, 0.75E13,& 0.10E13, 0.10E13, 0.10E13, 0.10E13 / ! use GEOS5/MERRA values for VGZ2, reichle+qliu, 8 Oct 2008 ! DATA VGZ2 /35.0, 20.0, 17.0, 0.6, 5.0, 0.6, 0.1, 0.1, 0.1, 0.1/ DATA VGZ2 /35.0, 20.0, 17.0, 0.6, 0.5, 0.6, 0.01, 0.01, 0.01, 0.01/ ! Dorman and Sellers (1989) DATA VGRT / 19700., 7000., 9400., 7000., 7000., 14000., 1., 1., 1., 1./ ! ----------------------------------------------------------------------------- DO ChNo = 1, NCH VEG = ITYP(ChNo) if(VEG.ne.9) then ! COMPUTE SOME PARAMETERS DIRECTLY: ! LAI and type dependent parameters; RDC formulation can use veg frac in next version. ! see the last line in the do loop - Sarith ! RDC (ChNo) = WLO * VGRDC (KLOW, ITYP (ChNo)) +& ! WHI * VGRDC (KLOW + 1, ITYP (ChNo)) ! use GEOS5/MERRA formulation for ROOTL, reichle+qliu, 8 Oct 2008 !ROOTL = WLO * VGROTL (KLOW, ITYP (ChNo)) +& ! WHI * VGROTL (KLOW + 1, ITYP (ChNo)) ROOTL = VGRT(VEG) ! the equations below are identical to those in AGCM tag "Ganymed-4_0" ! (GEOS_CatchGridComp.F90) ! - reichle, 5 Dec 2014 VROOT = ROOTL * VGROCA (VEG) DUM1 = ALOG (VROOT / (1. - VROOT)) DUM2 = 1. / (8. * MAPL_PI * ROOTL) ALPHAF = DUM2 * (VROOT - 3. -2. * DUM1) RSOIL1 (ChNo) = & VGRDRS (VEG) / (ROOTL * VGROTD (VEG)) RSOIL2 (ChNo) = ALPHAF / VGROTD (VEG) SCAT = GREEN (ChNo) * & (VGTR11 (VEG) + VGRF11 (VEG)) + & (1. - GREEN (ChNo)) * & (VGTR12 (VEG) + VGRF12 (VEG)) SQSCAT (ChNo) = SQRT (1. - SCAT) Z2 (ChNo) = VGZ2 (VEG) ! Modified time varied roughness scheme (Sarith 3/22/2016) ZVG(ChNo) = (Z2(ChNo)- 0.5*(Z2(ChNo)-MIN_VEG_HEIGHT)*exp(-ZLT(ChNo))) Z0 (ChNo) = Z0_BY_ZVEG * ZVG(ChNo) D (ChNo) = D0_BY_ZVEG * ZVG(ChNo) ! zero-plane displacement height ! Correction to RDC formulation -Randy Koster, 4/1/2011 ! bad: RDC = max(VGRDA(VEG)*min(VGRDB(VEG),LAI ),0.001) ! correct: RDC = max(VGRDA(VEG)*min(1., LAI/VGRDB(VEG)),0.001) !RDC(ChNo) = max(VGRDA(ITYP (ChNo))*min(VGRDB(ITYP (ChNo)), & ! ZLT (ChNo)),0.001) RDC(ChNo) = max( VGRDA(VEG)*min( 1., ZLT(ChNo)/VGRDB(VEG) ), 0.001) endif ENDDO ! RETURN END SUBROUTINE PMONTH ! ***************************************************************************** end module process_cat ! ========================== EOF ====================================