SURFEX v7.3
General documentation of Surfex
 All Classes Files Functions Variables Typedefs
/home/dasprezs/EXPORT_v7_3/src/SURFEX/mode_dstmblutl.F90
Go to the documentation of this file.
00001 ! ######spl
00002 MODULE MODE_DSTMBLUTL ! [mdl] Mobilization utilities
00003 !
00004 USE YOMHOOK   ,ONLY : LHOOK,   DR_HOOK
00005 USE PARKIND1  ,ONLY : JPRB
00006 !
00007 IMPLICIT NONE
00008 !
00009 CONTAINS
00010 !
00011 !----------------------------------------------------------------------------------------
00012 SUBROUTINE WND_FRC_THR_SLT_GET(PDNS_MDP, PDP, PWND_FRC_THR_SLT)  
00013 ! Purpose: Compute dry threshold friction velocity for saltation
00014 !++grini use dstgrd ! [mdl] Dust grid sizes
00015 !++grini use pmgrid ! [mdl] Spatial resolution parameters
00016 !++grini use dstcst ! [mdl] Physical constants for dust routines
00017 USE MODD_DSTMBL,    ONLY : XDNS_SLT
00018 USE MODD_CSTS,      ONLY : XG
00019 USE MODI_ABOR1_SFX
00020 !
00021 IMPLICIT NONE
00022 !
00023 ! Input
00024 !++grini real,intent(in)::dns_aer(dst_nbr) ! [kg m-3] Particle density
00025 !++grini real,intent(in)::dmt_aer(dst_nbr) ! [m] Particle diameter
00026 REAL, DIMENSION(:), INTENT(IN)   :: PDNS_MDP        ! [kg m-3] Midlayer density
00027 REAL, INTENT(IN)  :: PDP   
00028 ! Output
00029 REAL, DIMENSION(:), INTENT(OUT) :: PWND_FRC_THR_SLT ! [m s-1] Threshold friction velocity for saltation
00030 !
00031 ! Local
00032 REAL :: ZRYN, ZRYN1 ! [frc] Threshold friction Reynolds number approximation for optimal size
00033 REAL :: ZDNS_FCT    ! Density ratio factor for saltation calculation
00034 REAL :: ZICF_FCT    ! Interparticle cohesive forces factor for saltation calculation
00035 REAL :: ZTMP        ! Factor in saltation computation
00036 INTEGER :: I        ! [idx] Counting index
00037 !
00038 REAL(KIND=JPRB) :: ZHOOK_HANDLE   
00039 !
00040 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:WND_FRC_THR_SLT_GET',0,ZHOOK_HANDLE)
00041 !
00042 ! Initialize some variables
00043 ! MaB95 pzn. for Re*t(D_opt) circumvents iterative solution
00044 ZRYN = 0.38d0 + 1331.0d0*(100.0d0*PDP)**1.56d0                        ! [frc] "B" MaB95 p. 16417 (5)
00045 ! Given Re*t(D_opt), compute time independent factors contributing to u*t
00046 ZDNS_FCT                 = XDNS_SLT * XG * PDP                        ! IvW82 p. 115 (6) MaB95 p. 16417 (4)
00047 ZICF_FCT                 = 1.0d0 + 6.0d-07/(XDNS_SLT*XG*(PDP**2.5d0)) ! [frc] IvW82 p. 115 (6) MaB95 p. 16417 (4) Interparticle cohesive forces
00048 !
00049 IF (ZRYN < 0.03d0) THEN
00050   CALL ABOR1_SFX('MODE_DSTMBLUTL:WND_FRC_THR_SLT_GET: RYN < 0.03')
00051 ELSEIF (ZRYN < 10.0d0) THEN
00052   ZRYN1 = -1.0d0 + 1.928d0 * (ZRYN**0.0922d0) ! [frc] IvW82 p. 114 (3), MaB95 p. 16417 (6)
00053   ZRYN1 = 0.1291d0 * 0.1291d0 / ZRYN1          ! [frc] 
00054 ELSE
00055   ZRYN1 = 1.0d0 - 0.0858d0 * EXP(-0.0617d0*(ZRYN-10.0d0)) ! [frc] IvW82 p. 114 (3), MaB95 p. 16417 (7)
00056   ZRYN1 = 0.120d0 **2 * ZRYN1**2               ! [frc]
00057   !ryn1=0.129*0.129*ryn1*ryn1 dans le cas mm, à vérifier
00058 ENDIF 
00059 !
00060 ! This method minimizes the number of square root computations performed
00061 ZTMP = SQRT (ZICF_FCT*ZDNS_FCT*ZRYN1)
00062 DO I = 1, SIZE(PDNS_MDP)
00063   PWND_FRC_THR_SLT(I) = ZTMP / SQRT(PDNS_MDP(I))
00064 ENDDO
00065 !
00066 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:WND_FRC_THR_SLT_GET',1,ZHOOK_HANDLE)
00067 !
00068 END SUBROUTINE WND_FRC_THR_SLT_GET
00069 !
00070 !----------------------------------------------------------------------------------------
00071 SUBROUTINE VWC2GWC (OFLG_MBL, PVWC_SAT, PVWC_SFC, PGWC_SFC)
00072 ! Purpose: Convert volumetric water content to gravimetric water content
00073 !++alfgr use pmgrid ! [mdl] Spatial resolution parameters
00074 use MODD_CSTS, only : XRHOLI  ! [kg/m3] density of liquid water
00075 USE MODD_DSTMBL, ONLY : XDNS_SLT
00076 !
00077 !++alfgr use dstblm,only:dns_H2O_lqd_std ! [mdl] Boundary layer meteorology for non-vegetated land surfaces
00078 IMPLICIT NONE
00079 ! Input
00080 LOGICAL, DIMENSION(:), INTENT(IN) :: OFLG_MBL ! [flg] Mobilization candidate flag
00081 REAL,    DIMENSION(:), INTENT(IN) :: PVWC_SAT ! [m3 m-3] Saturated volumetric water content (sand-dependent)
00082 REAL,    DIMENSION(:), INTENT(IN) :: PVWC_SFC ! [m3 m-3] Volumetric water content
00083 REAL,    DIMENSION(:), INTENT(OUT):: PGWC_SFC ! [kg kg-1] Gravimetric water content
00084 ! Local
00085 REAL, DIMENSION(SIZE(PVWC_SAT)) :: ZDNS_BLK_DRY ! [kg m-3] Bulk density of dry surface soil
00086 INTEGER :: I
00087 REAL(KIND=JPRB) :: ZHOOK_HANDLE 
00088 ! Main Code
00089 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:VWC2GWC',0,ZHOOK_HANDLE)
00090 !
00091 ! Initialize output
00092 DO I=1,SIZE(PVWC_SAT)
00093   IF (OFLG_MBL(I)) THEN
00094     ! Assume volume of air pores when dry equals saturated VWC
00095     ! This implies air pores are completely filled by water in saturated soil
00096     ZDNS_BLK_DRY(I) = XDNS_SLT * (1.0 - PVWC_SAT(I)) ! [kg m-3] Bulk density of dry surface soil
00097     PGWC_SFC    (I) = PVWC_SFC(I) * XRHOLI / ZDNS_BLK_DRY(I)
00098   ENDIF
00099 ENDDO
00100 !
00101 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:VWC2GWC',1,ZHOOK_HANDLE)
00102 END SUBROUTINE VWC2GWC
00103 !
00104 !----------------------------------------------------------------------------------------  
00105 SUBROUTINE FRC_THR_NCR_WTR_GET(OFLG_MBL, PGWC_THR, PGWC_SFC, PFRC_THR_NCR_WTR)   
00106 ! Purpose: Compute factor by which soil moisture increases threshold friction velocity 
00107 ! This parameterization is based on FMB99
00108 !++alfgr use pmgrid ! [mdl] Spatial resolution parameters
00109 IMPLICIT NONE
00110 ! dans le calcul de l'humidit? du sol(a posteriori)
00111 ! Input
00112 LOGICAL, DIMENSION(:), INTENT(IN) :: OFLG_MBL         ! I [flg] Mobilization candidate flag
00113 REAL,    DIMENSION(:), INTENT(IN) :: PGWC_THR         ! [kg kg-1] Threshold gravimetric water content
00114 REAL,    DIMENSION(:), INTENT(IN) :: PGWC_SFC         ! [kg kg-1] Gravimetric water content
00115 REAL,    DIMENSION(:), INTENT(OUT):: PFRC_THR_NCR_WTR ! [frc] Factor by which moisture increases threshold friction velocity
00116 ! Local
00117 INTEGER :: I ! [idx] Counting index
00118 REAL(KIND=JPRB) :: ZHOOK_HANDLE 
00119 
00120 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FRC_THR_NCR_WTR_GET',0,ZHOOK_HANDLE)
00121 ! Main Code
00122 ! Initialize output
00123 PFRC_THR_NCR_WTR(:) = 1.0d0 ! [frc] Factor by which moisture increases threshold friction velocity
00124 !
00125 DO I = 1, SIZE(PGWC_SFC)
00126   IF (OFLG_MBL(I)) THEN
00127     ! Adjust threshold velocity for inhibition by moisture
00128     ! frc_thr_ncr_wtr(lon_idx)=exp(22.7d0*vwc_sfc(lon_idx)) ! [frc] SRL96
00129     ! Compute threshold soil moisture based on clay content
00130     IF (PGWC_SFC(I) > PGWC_THR(I)) &
00131       PFRC_THR_NCR_WTR(I) = SQRT(1.0d0 + 1.21d0 * (100.0d0*(PGWC_SFC(I)-PGWC_THR(I)))**0.68d0) ! [frc] FMB99 p. 155 (15)
00132   ENDIF
00133 ENDDO
00134 ! Uncomment following line to remove all dependence on gwc_sfc
00135 ! frc_thr_ncr_wtr(lon_idx)=1.0d0 ! [frc]
00136 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FRC_THR_NCR_WTR_GET',1,ZHOOK_HANDLE)
00137 END SUBROUTINE FRC_THR_NCR_WTR_GET
00138 !  
00139 !----------------------------------------------------------------------------------------  
00140 !++alfgr fxm: Fix this so that we can actually use rgh_mmn_mbl different for grid cells
00141 SUBROUTINE FRC_THR_NCR_DRG_GET(PRGH_MMN_MBL, PRGH_MMN_SMT, PFRC_THR_NCR_DRG)
00142 ! Purpose: Compute factor by which surface roughness increases threshold friction velocity
00143 ! This parameterization is based on MaB95 and GMB98
00144 !++alfgr use pmgrid ! [mdl] Spatial resolution parameters
00145 USE MODI_ABOR1_SFX
00146 !
00147 IMPLICIT NONE
00148 ! Input
00149 REAL, INTENT(IN) :: PRGH_MMN_MBL     ! [m] Roughness length momentum for erodible surfaces
00150 REAL, INTENT(IN) :: PRGH_MMN_SMT                   ! [m] Smooth roughness length
00151 ! Output
00152 REAL, INTENT(OUT):: PFRC_THR_NCR_DRG ! [frc] Factor by which roughness increases threshold friction velocity
00153 ! Local
00154 REAL :: ZWND_FRC_FSH_FRC     ! [frc] Efficient fraction of wind friction
00155 REAL(KIND=JPRB) :: ZHOOK_HANDLE 
00156 !
00157 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FRC_THR_NCR_DRG_GET',0,ZHOOK_HANDLE)
00158 ! Main Code
00159 !
00160 ! Adjust threshold velocity for inhibition by roughness elements
00161 ZWND_FRC_FSH_FRC = & ! [frc] MaB95 p. 16420, GMB98 p. 6207
00162   1.0d0 - LOG(PRGH_MMN_MBL/PRGH_MMN_SMT) / LOG(0.35d0 * ((0.1d0/PRGH_MMN_SMT)**0.8d0))
00163 ZWND_FRC_FSH_FRC = MAX(1.e-6, MIN(1., ZWND_FRC_FSH_FRC))
00164 IF (ZWND_FRC_FSH_FRC <= 0.0d0 .OR. ZWND_FRC_FSH_FRC > 1.0d0) &
00165   CALL ABOR1_SFX("MODE_DSTMBLUTL:FRC_THR_NCR_DRG_GET0: WND_FRC_FSH_FRC OUT OF RANGE")
00166 PFRC_THR_NCR_DRG = 1.0d0 / ZWND_FRC_FSH_FRC! [frc]
00167 ! fxm: 19991012 Set frc_thr_ncr_drg=1.0, equivalent to assuming mobilization takes place at smooth roughness length
00168 !++alfgr removed this frc_thr_ncr_drg=1.0d0
00169 !
00170 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FRC_THR_NCR_DRG_GET',1,ZHOOK_HANDLE)
00171 END SUBROUTINE FRC_THR_NCR_DRG_GET
00172 !
00173 !----------------------------------------------------------------------------------------
00174 SUBROUTINE WND_FRC_SLT_GET(OFLG_MBL, PWND_FRC, PWND_RFR, PWND_RFR_THR_SLT, PWND_FRC_SLT)
00175 ! Purpose: Compute the saltating friction velocity
00176 ! Saltation increases friction speed by roughening surface, AKA "Owen's effect"
00177 ! This acts as a positive feedback to the friction speed
00178 ! GMB98 parameterized this feedback in terms of 10 m windspeeds
00179 !++alfgr use pmgrid ! [mdl] Spatial resolution parameters
00180 IMPLICIT NONE
00181 ! Input
00182 LOGICAL, DIMENSION(:), INTENT(IN) :: OFLG_MBL         ! I [flg] Mobilization candidate flag
00183 REAL,    DIMENSION(:), INTENT(IN) :: PWND_FRC         ! I [m s-1] Surface friction velocity
00184 REAL,    DIMENSION(:), INTENT(IN) :: PWND_RFR         ! I [m s-1] Wind speed at reference height
00185 REAL,    DIMENSION(:), INTENT(IN) :: PWND_RFR_THR_SLT ! I [m s-1] Threshold 10 m wind speed for saltation
00186 REAL,    DIMENSION(:), INTENT(OUT):: PWND_FRC_SLT     ! O [m s-1] Saltating friction velocity
00187 ! Local
00188 REAL :: ZWND_RFR_DLT     ! [m s-1] Reference windspeed excess over threshold
00189 REAL :: ZWND_FRC_SLT_DLT ! [m s-1] Friction velocity increase from saltation
00190 INTEGER :: I             ! [idx] Counting index
00191 REAL(KIND=JPRB) :: ZHOOK_HANDLE 
00192 !
00193 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:WND_FRC_SLT_GET',0,ZHOOK_HANDLE)   
00194 ! Main Code
00195 ! Compute saltating friction velocity, accounting for "Owen's effect"
00196 PWND_FRC_SLT(:) = PWND_FRC(:) ! [m s-1] Saltating friction velocity
00197 !
00198 DO I = 1, SIZE(PWND_FRC)
00199   IF (OFLG_MBL(I) .AND. PWND_RFR(I) >= PWND_RFR_THR_SLT(I)) THEN
00200     ! Saltation roughens the boundary layer, AKA "Owen's effect"
00201     ! GMB98 p. 6206 Fig. 1 shows observed/computed u* dependence on observed U(1 m)
00202     ! GMB98 p. 6209 (12) has u* in cm s-1 and U, Ut in m s-1, personal communication, D. Gillette, 19990529
00203     ! With everything in MKS, the 0.3 coefficient in GMB98 (12) becomes 0.003 
00204     ! Increase in friction velocity due to saltation varies as square of 
00205     ! difference between reference wind speed and reference threshold speed 
00206     ZWND_RFR_DLT     = PWND_RFR(I) - PWND_RFR_THR_SLT(I)
00207     ZWND_FRC_SLT_DLT = 0.003d0 * ZWND_RFR_DLT**2          ! [m s-1] Friction velocity increase from saltation GMB98 p. 6209
00208     PWND_FRC_SLT(I)  = PWND_FRC_SLT(I) + ZWND_FRC_SLT_DLT ! [m s-1] Saltating friction velocity
00209   ENDIF
00210 ENDDO
00211 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:WND_FRC_SLT_GET',1,ZHOOK_HANDLE)
00212 END SUBROUTINE WND_FRC_SLT_GET
00213 !
00214 !----------------------------------------------------------------------------------------
00215 SUBROUTINE FLX_MSS_HRZ_SLT_TTL_WHI79_GET(PCOEFF, OFLG_MBL, PDNS_MDP, PWND_FRC, &
00216                         PWND_FRC_THR_SLT, PFLX_MSS_HRZ_SLT_TTL)
00217 ! Purpose: Compute vertically integrated streamwise mass flux of particles
00218 ! Theory: Uses method proposed by White (1979)
00219 ! fxm: use surface air density not midlayer density
00220 !++alfgr use pmgrid ! [mdl] Spatial resolution parameters
00221 !++alfgr use dstcst ! [mdl] Physical constants for dust routines
00222 use MODD_CSTS,      ONLY : XG  ! Gravitation constant
00223 !
00224 IMPLICIT NONE
00225 ! Input
00226 REAL,    DIMENSION(:), INTENT(IN) :: PCOEFF
00227 LOGICAL, DIMENSION(:), INTENT(IN) :: OFLG_MBL             ! I [flg] Mobilization candidate flag
00228 REAL,    DIMENSION(:), INTENT(IN) :: PDNS_MDP             ! I [kg m-3] Midlayer density
00229 REAL,    DIMENSION(:), INTENT(IN) :: PWND_FRC             ! I [m s-1] Surface friction velocity
00230 REAL,    DIMENSION(:), INTENT(IN) :: PWND_FRC_THR_SLT     ! I [m s-1] Threshold friction speed for saltation
00231 REAL,    DIMENSION(:), INTENT(OUT):: PFLX_MSS_HRZ_SLT_TTL ! O [kg m-1 s-1] Vertically integrated streamwise mass flux
00232 ! Local
00233 REAL :: ZWND_FRC_RAT ! [frc] Ratio of wind friction threshold to wind friction
00234 INTEGER :: I    ! [idx] Counting index for lon
00235 REAL(KIND=JPRB) :: ZHOOK_HANDLE 
00236 !
00237 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FLX_MSS_HRZ_SLT_TTL_WHI79_GET',0,ZHOOK_HANDLE)   
00238 ! Main Code
00239 ! Initialize output
00240 PFLX_MSS_HRZ_SLT_TTL(:) = 0.0d0 ! [kg m-1 s-1]
00241 !
00242 DO I = 1, SIZE(PDNS_MDP)
00243   IF (OFLG_MBL(I) .AND. PWND_FRC(I) > PWND_FRC_THR_SLT(I)) THEN
00244     ZWND_FRC_RAT = PWND_FRC_THR_SLT(I) / PWND_FRC(I)            ! [frc]
00245     PFLX_MSS_HRZ_SLT_TTL(I) =                                 & ! [kg m-1 s-1]
00246       PCOEFF(I) * PDNS_MDP(I) * (PWND_FRC(I)**3.0d0) *      &
00247       (1.0d0 - ZWND_FRC_RAT) * (1.0d0 + ZWND_FRC_RAT)**2 / XG   ! Whi79 p. 4648 (19), MaB97 p. 16422 (28)
00248   ENDIF
00249 ENDDO
00250 
00251 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FLX_MSS_HRZ_SLT_TTL_WHI79_GET',1,ZHOOK_HANDLE) 
00252 END SUBROUTINE FLX_MSS_HRZ_SLT_TTL_WHI79_GET
00253 ! 
00254 !----------------------------------------------------------------------------------------
00255 SUBROUTINE FLX_MSS_VRT_DST_TTL_MAB95_GET(OFLG_MBL, PMSS_FRC_CLY, PFLX_MSS_HRZ_SLT_TTL, & 
00256                         PDST_SLT_FLX_RAT_TTL, PFLX_MSS_VRT_DST_TTL)
00257 ! Purpose: Diagnose total vertical mass flux of dust from vertically integrated streamwise mass flux
00258 ! Theory: Uses clay-based method proposed by Marticorena & Bergametti (1995)
00259 ! Their parameterization is based only on data for mss_frc_cly < 0.20
00260 ! For clayier soils, dst_slt_flx_rat_ttl may behave dramatically differently
00261 ! Whether this behavior changes when mss_frc_cly > 0.20 is unknown
00262 ! Anecdotal evidence suggests vertical flux decreases for mss_frc_cly > 0.20
00263 ! Thus we use min[mss_frc_cly,0.20] in MaB95 parameterization
00264 IMPLICIT NONE
00265 ! Input
00266 LOGICAL, DIMENSION(:), INTENT(IN) :: OFLG_MBL             ! I [flg] Mobilization candidate flag
00267 REAL,    DIMENSION(:), INTENT(IN) :: PMSS_FRC_CLY         ! I [frc] Mass fraction clay
00268 REAL,    DIMENSION(:), INTENT(IN) :: PFLX_MSS_HRZ_SLT_TTL ! I [kg m-1 s-1] Vertically integrated streamwise mass flux
00269 REAL,    DIMENSION(:), INTENT(OUT):: PDST_SLT_FLX_RAT_TTL ! O [m-1] Ratio of vertical dust flux to streamwise mass flux
00270 REAL,    DIMENSION(:), INTENT(OUT):: PFLX_MSS_VRT_DST_TTL ! O [kg m-2 s-1] Total vertical mass flux of dust
00271 ! Local
00272 REAL :: ZMSS_FRC_CLY_VLD ! [frc] Mass fraction clay limited to 0.20
00273 INTEGER :: I             ! [idx] Counting index for lon
00274 REAL(KIND=JPRB) :: ZHOOK_HANDLE 
00275 !
00276 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FLX_MSS_VRT_DST_TTL_MAB95_GET',0,ZHOOK_HANDLE)  
00277 !
00278 DO I = 1, SIZE(PMSS_FRC_CLY)
00279   IF (OFLG_MBL(I)) THEN
00280     ! 19990603: fxm: Dust production is EXTREMELY sensitive to this parameter, which changes flux by 3 orders of magnitude in 0.0 < mss_frc_cly < 0.20
00281     ZMSS_FRC_CLY_VLD = MIN(PMSS_FRC_CLY(I), 0.2) ! [frc]
00282     PDST_SLT_FLX_RAT_TTL(I) = & ! [m-1]
00283       100.0d0 * EXP(LOG(10.0d0)*(13.4d0*ZMSS_FRC_CLY_VLD - 6.0d0))              ! MaB95 p. 16423 (47)
00284     PFLX_MSS_VRT_DST_TTL(I) = PFLX_MSS_HRZ_SLT_TTL(I) * PDST_SLT_FLX_RAT_TTL(I) ! [kg m-1 s-1]
00285   ENDIF
00286 ENDDO
00287 !
00288 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FLX_MSS_VRT_DST_TTL_MAB95_GET',1,ZHOOK_HANDLE) 
00289 END SUBROUTINE FLX_MSS_VRT_DST_TTL_MAB95_GET
00290 !
00291 !----------------------------------------------------------------------------------------
00292 SUBROUTINE FLX_MSS_VRT_DST_AUST_GET(OFLG_MBL, PDNS_MDP, PWND_FRC_THR_SLT, &
00293         PFLX_MSS_HRZ_SLT_TTL, PDST_SLT_FLX_RAT_TTL, PFLX_MSS_VRT_DST_TTL)
00294 ! Purpose: Diagnose total vertical mass flux of dust from vertically integrated streamwise mass flux
00295 use MODD_CSTS,      ONLY : XG  ! Gravitation constant
00296 USE MODD_DSTMBL, ONLY : XDMT_SLT_OPT, XDMT_ERO_OPT, XDNS_SLT, XGAMA
00297 !
00298 IMPLICIT NONE
00299 ! Input
00300 LOGICAL, DIMENSION(:), INTENT(IN) :: OFLG_MBL             ! I [flg] Mobilization candidate flag
00301 REAL,    DIMENSION(:), INTENT(IN) :: PDNS_MDP             ! I [kg m-3] Midlayer density
00302 REAL,    DIMENSION(:), INTENT(IN) :: PWND_FRC_THR_SLT     ! I [m s-1] Threshold friction speed for saltation
00303 REAL,    DIMENSION(:), INTENT(IN) :: PFLX_MSS_HRZ_SLT_TTL ! I [kg m-1 s-1] Vertically integrated streamwise mass flux
00304 REAL,    DIMENSION(:), INTENT(OUT):: PDST_SLT_FLX_RAT_TTL ! O [m-1] Ratio of vertical dust flux to streamwise mass flux
00305 REAL,    DIMENSION(:), INTENT(OUT):: PFLX_MSS_VRT_DST_TTL ! O [kg m-2 s-1] Total vertical mass flux of dust
00306 ! Local
00307 REAL, DIMENSION(SIZE(PDNS_MDP)) :: ZROP_ROA  ! [frc] ratio of particle densite to Midlayer density
00308 REAL :: ZEXPDD, ZLNDS, ZBETA, ZBGXG
00309 INTEGER :: I
00310 REAL(KIND=JPRB) :: ZHOOK_HANDLE 
00311 !
00312 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FLX_MSS_VRT_DST_AUST_GET',0,ZHOOK_HANDLE)  
00313 !
00314 ! Initialize some variables
00315 ZEXPDD = EXP(-140.7d0 * XDMT_ERO_OPT + 0.37d0)
00316 ZLNDS  = (0.328*1.0d-4) + (0.125d0*1.0d-4*LOG(XDMT_SLT_OPT*1.0E+3))
00317 ZBETA  = ZLNDS * ZEXPDD ! [mm]
00318 ZBGXG  = ZBETA * XGAMA * XG
00319 !
00320 DO I = 1, SIZE(PDNS_MDP)
00321   IF (OFLG_MBL(I)) THEN
00322     ZROP_ROA(I) = XDNS_SLT / PDNS_MDP(I) 
00323     PDST_SLT_FLX_RAT_TTL(I) =          & ! [mm/m]
00324       (2.0d0/3.0d0) * ZBGXG * ZROP_ROA(I) / (PWND_FRC_THR_SLT(I)**2.0d0)
00325     PFLX_MSS_VRT_DST_TTL(I) = 1.0d-3 * & ! [kg m-1 s-1]
00326       PFLX_MSS_HRZ_SLT_TTL(I) * PDST_SLT_FLX_RAT_TTL(I)
00327   ENDIF
00328 ENDDO
00329 !
00330 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:FLX_MSS_VRT_DST_AUST_GET',1,ZHOOK_HANDLE)
00331 END SUBROUTINE FLX_MSS_VRT_DST_AUST_GET
00332 !
00333 !----------------------------------------------------------------------------------------
00334 SUBROUTINE DISTRIBUTION(KMODE, PDLNDP, KTEXT, PDP, PDSRLV, PZS0)
00335 !
00336 USE MODD_CSTS,      ONLY : XPI
00337 USE MODD_DSTMBL, ONLY : XDNS_SLT
00338 !
00339 IMPLICIT NONE
00340 !
00341 !INPUT, set their dimensions to their passed lengths or to KSIZE ?
00342 INTEGER, INTENT(IN) :: KMODE  ![nbr] loop over mode
00343 REAL,    INTENT(IN) :: PDLNDP ! delta Dp [?m]
00344 INTEGER, DIMENSION(:), INTENT(IN) :: KTEXT  ! texture
00345 REAL, DIMENSION(:),   INTENT(OUT) :: PDP
00346 REAL, DIMENSION(:,:), INTENT(OUT) :: PDSRLV ! [--] Output surface relative
00347 REAL, DIMENSION(:),   INTENT(OUT) :: PZS0   ! [m] Output rugosit? de la surface lisse (Dp/30)
00348 !Define local variables:
00349 REAL, DIMENSION(KMODE,12) :: ZSIGMA ! déviation géométrique standard du mode
00350 REAL, DIMENSION(KMODE,12) :: ZMFRAC ! fraction massique du mode
00351 REAL, DIMENSION(KMODE,12) :: ZDMED  ! diametre median du mode
00352 REAL, DIMENSION(SIZE(PDSRLV,1),SIZE(PDSRLV,2)) :: ZDS
00353 REAL, DIMENSION(SIZE(PDSRLV,2)) :: ZDPLN, ZSP, ZDMLN
00354 REAL :: ZDM1, ZDM2
00355 INTEGER :: IMOD  ! Counter for number of mode
00356 INTEGER :: IDP   ! Counter for number of particle
00357 INTEGER :: ITEX  ! Counter for number of texture
00358 REAL(KIND=JPRB) :: ZHOOK_HANDLE 
00359 !
00360 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:DISTRIBUTION',0,ZHOOK_HANDLE)  
00361 !
00362 !1 Sand
00363 !2 Loamy Sand
00364 !3 Sandy Loam
00365 !4 Silt Loam
00366 !5 Loam
00367 !6 Sandy Clay Loam
00368 !7 Silty Clay Loam
00369 !8 Clay Loam
00370 !9 Sandy Clay
00371 !10 Silty Clay
00372 !11 Clay
00373 !12 Silt
00374 ZMFRAC(1,:) = (/0.0, 0.1, 0.1, 0.15, 0.15, 0.2, 0.2, 0.3, 0.35, 0.4, 0.5, 0.15/)
00375 ZMFRAC(2,:) = (/0.1, 0.3, 0.3, 0.35, 0.50, 0.5, 0.5, 0.5, 0.00, 0.0, 0.0, 0.40/)
00376 ZMFRAC(3,:) = (/0.9, 0.6, 0.6, 0.50, 0.35, 0.3, 0.3, 0.2, 0.65, 0.6, 0.5, 0.45/)
00377 !
00378 ZSIGMA(1,:) = (/1.8, 1.8, 1.8, 1.8, 1.8, 1.8, 1.8, 1.8, 1.8, 1.8, 1.8, 1.8/)
00379 ZSIGMA(2,:) = (/1.7, 1.7, 1.7, 1.7, 1.7, 1.7, 1.7, 1.7, 1.8, 1.8, 1.8, 1.7/)
00380 ZSIGMA(3,:) = (/1.6, 1.6, 1.6, 1.6, 1.6, 1.7, 1.7, 1.7, 1.8, 1.8, 1.8, 1.6/)
00381 !
00382 ZDMED (1,:) = (/  10., 10.,  5.,  5., 2.5, 2.5, 2.5,  1.,  1., 0.5, 0.5, 2.5/)
00383 ZDMED (2,:) = (/ 100.,100.,100.,100., 75., 75., 50., 50., 10., 10., 10., 75./)
00384 ZDMED (3,:) = (/1000.,690.,520.,520.,520.,210.,210.,125.,100.,100.,100.,520./)
00385 !
00386 DO ITEX = 1, SIZE(PZS0) 
00387   PZS0(ITEX) =1E-6*ZDMED(3,ITEX)/30.0
00388 ENDDO
00389 !
00390 ! Initialisation pour IDP = 1
00391 PDP  (1) = 0.1          ! ?m
00392 ZDPLN(1) = LOG(PDP(1))  ! ?m
00393 ! 
00394 DO IDP = 2, SIZE(PDSRLV,2)
00395   ZDPLN(IDP) = ZDPLN(IDP-1) + PDLNDP
00396   PDP  (IDP) = EXP(ZDPLN(IDP))
00397 END DO
00398 !
00399 ! calcul pour chaque particule 
00400 ! Boucle NDP 
00401 ! Initialisation
00402 ! surface basale totale
00403 ZSP(:) = 0.d0
00404 ! Calcul de la distribution massique des particules (dm(Dp)/dln(Dp)) pour chaque particule de diametre Dp
00405 DO ITEX = 1, SIZE(PDSRLV,1)
00406   DO IDP = 1, SIZE(PDSRLV,2)
00407     ZDMLN(IDP) = 0.0d0
00408     DO IMOD = 1, KMODE
00409       ZDM1 = ZMFRAC(IMOD,ITEX) / (SQRT(2.d0*XPI) * LOG(ZSIGMA(IMOD,ITEX)))
00410       ZDM2 = EXP((LOG(PDP(IDP)) - LOG(ZDMED(IMOD,ITEX)))**2. / (-2.d0 * (LOG(ZSIGMA(IMOD,ITEX)))**2.d0))
00411       ZDMLN(IDP) = ZDMLN(IDP) + ZDM1*ZDM2
00412     END DO 
00413     ! Calcul de la surface basale occup?e par chaque particule de diam?tre Dp  
00414     ZDS(ITEX,IDP) = 3.d0 * ZDMLN(IDP) * PDLNDP / (2.* XDNS_SLT * PDP(IDP))
00415     ! Calcul de la surface basale totale Sp= (somme (ds(Dp)):
00416     ZSP(ITEX) = ZSP(ITEX) + ZDS(ITEX,IDP)
00417   ENDDO
00418 ENDDO 
00419 !
00420 ! Calcul de la distribution de surface relative des particules
00421 DO ITEX = 1, SIZE(PZS0)
00422   DO IDP = 1, SIZE(PDSRLV,2)
00423     PDSRLV(ITEX,IDP) = ZDS(ITEX,IDP) / ZSP(ITEX)
00424     ENDDO 
00425 END DO 
00426 ! 
00427 IF (LHOOK) CALL DR_HOOK('MODE_DSTMBLUTL:DISTRIBUTION',1,ZHOOK_HANDLE)  
00428 END SUBROUTINE DISTRIBUTION
00429 !
00430 END MODULE MODE_DSTMBLUTL