!WRF:MEDIATION_LAYER:FIRE_MODEL
! Routines dealing with the atmosphere
module module_fr_fire_atm 1
use module_model_constants
, only: cp,xlv
use module_fr_fire_util
contains
SUBROUTINE fire_tendency( & 1,4
ids,ide, kds,kde, jds,jde, & ! dimensions
ims,ime, kms,kme, jms,jme, &
its,ite, kts,kte, jts,jte, &
grnhfx,grnqfx,canhfx,canqfx, & ! heat fluxes summed up to atm grid
alfg,alfc,z1can, & ! coeffients, properties, geometry
zs,z_at_w,dz8w,mu,rho, &
rthfrten,rqvfrten) ! theta and Qv tendencies
! This routine is atmospheric physics
! it does NOT go into module_fr_fire_phys because it is not related to fire physical processes
! --- this routine takes fire generated heat and moisture fluxes and
! calculates their influence on the theta and water vapor
! --- note that these tendencies are valid at the Arakawa-A location
IMPLICIT NONE
! --- incoming variables
INTEGER , INTENT(in) :: ids,ide, kds,kde, jds,jde, &
ims,ime, kms,kme, jms,jme, &
its,ite, kts,kte, jts,jte
REAL, INTENT(in), DIMENSION( ims:ime,jms:jme ) :: grnhfx,grnqfx ! W/m^2
REAL, INTENT(in), DIMENSION( ims:ime,jms:jme ) :: canhfx,canqfx ! W/m^2
REAL, INTENT(in), DIMENSION( ims:ime,jms:jme ) :: zs ! topography (m abv sealvl)
REAL, INTENT(in), DIMENSION( ims:ime,jms:jme ) :: mu ! dry air mass (Pa)
REAL, INTENT(in), DIMENSION( ims:ime,kms:kme,jms:jme ) :: z_at_w ! m abv sealvl
REAL, INTENT(in), DIMENSION( ims:ime,kms:kme,jms:jme ) :: dz8w ! dz across w-lvl
REAL, INTENT(in), DIMENSION( ims:ime,kms:kme,jms:jme ) :: rho ! density
REAL, INTENT(in) :: alfg ! extinction depth surface fire heat (m)
REAL, INTENT(in) :: alfc ! extinction depth crown fire heat (m)
REAL, INTENT(in) :: z1can ! height of crown fire heat release (m)
! --- outgoing variables
REAL, INTENT(out), DIMENSION( ims:ime,kms:kme,jms:jme ) :: &
rthfrten, & ! theta tendency from fire (in mass units)
rqvfrten ! Qv tendency from fire (in mass units)
! --- local variables
INTEGER :: i,j,k
INTEGER :: i_st,i_en, j_st,j_en, k_st,k_en
REAL :: cp_i
REAL :: rho_i
REAL :: xlv_i
REAL :: z_w
REAL :: fact_g, fact_c
REAL :: alfg_i, alfc_i
REAL, DIMENSION( its:ite,kts:kte,jts:jte ) :: hfx,qfx
!! character(len=128)::msg
do j=jts,jte
do k=kts,min(kte+1,kde)
do i=its,ite
rthfrten(i,k,j)=0.
rqvfrten(i,k,j)=0.
enddo
enddo
enddo
! --- set some local constants
cp_i = 1./cp ! inverse of specific heat
xlv_i = 1./xlv ! inverse of latent heat
alfg_i = 1./alfg
alfc_i = 1./alfc
!!write(msg,'(8e11.3)')cp,cp_i,xlv,xlv_i,alfg,alfc,z1can
!!call message(msg)
call print_2d_stats
(its,ite,jts,jte,ims,ime,jms,jme,grnhfx,'fire_tendency:grnhfx')
call print_2d_stats
(its,ite,jts,jte,ims,ime,jms,jme,grnqfx,'fire_tendency:grnqfx')
! --- set loop indicies : note that
i_st = MAX(its,ids+1)
i_en = MIN(ite,ide-1)
k_st = kts
k_en = MIN(kte,kde-1)
j_st = MAX(jts,jds+1)
j_en = MIN(jte,jde-1)
! --- distribute fluxes
DO j = j_st,j_en
DO k = k_st,k_en
DO i = i_st,i_en
! --- set z (in meters above ground)
z_w = z_at_w(i,k,j) - zs(i,j) ! should be zero when k=k_st
! --- heat flux
fact_g = cp_i * EXP( - alfg_i * z_w )
IF ( z_w < z1can ) THEN
fact_c = cp_i
ELSE
fact_c = cp_i * EXP( - alfc_i * (z_w - z1can) )
END IF
hfx(i,k,j) = fact_g * grnhfx(i,j) + fact_c * canhfx(i,j)
!! write(msg,2)i,k,j,z_w,grnhfx(i,j),hfx(i,k,j)
!!2 format('hfx:',3i4,6e11.3)
!! call message(msg)
! --- vapor flux
fact_g = xlv_i * EXP( - alfg_i * z_w )
IF (z_w < z1can) THEN
fact_c = xlv_i
ELSE
fact_c = xlv_i * EXP( - alfc_i * (z_w - z1can) )
END IF
qfx(i,k,j) = fact_g * grnqfx(i,j) + fact_c * canqfx(i,j)
!! if(hfx(i,k,j).ne.0. .or. qfx(i,k,j) .ne. 0.)then
!! write(msg,1)i,k,j,hfx(i,k,j),qfx(i,k,j)
!!1 format('tend:',3i6,2e11.3)
!! call message(msg)
! endif
END DO
END DO
END DO
! --- add flux divergence to tendencies
!
! multiply by dry air mass (mu) to eliminate the need to
! call sr. calculate_phy_tend (in dyn_em/module_em.F)
DO j = j_st,j_en
DO k = k_st,k_en-1
DO i = i_st,i_en
rho_i = 1./rho(i,k,j)
rthfrten(i,k,j) = - mu(i,j) * rho_i * (hfx(i,k+1,j)-hfx(i,k,j)) / dz8w(i,k,j)
rqvfrten(i,k,j) = - mu(i,j) * rho_i * (qfx(i,k+1,j)-qfx(i,k,j)) / dz8w(i,k,j)
END DO
END DO
END DO
call print_3d_stats
(its,ite,kts,kte,jts,jte,ims,ime,kms,kme,jms,jme,rthfrten,'fire_tendency:rthfrten')
call print_3d_stats
(its,ite,kts,kte,jts,jte,ims,ime,kms,kme,jms,jme,rqvfrten,'fire_tendency:rqvfrten')
RETURN
END SUBROUTINE fire_tendency
!
!***
!
end module module_fr_fire_atm