Skip to content

Commit

Permalink
split normal and advection motion of the fireline in tend_ls
Browse files Browse the repository at this point in the history
  • Loading branch information
Jan Mandel committed Sep 26, 2008
1 parent 86dfa31 commit 4b5d8b2
Show file tree
Hide file tree
Showing 2 changed files with 52 additions and 22 deletions.
72 changes: 51 additions & 21 deletions wrfv2_fire/phys/module_fr_sfire_core.F
Expand Up @@ -1143,7 +1143,7 @@ subroutine prop_ls( id, & ! for debug
IMTS,IMTE,JMTS,JMTE, & ! memory dims for tend
ids,ide,jds,jde, & ! domain dims - where lfn exists
ihs,ihe,jhs,jhe, & ! where tend computed
ts,dx,dy, & ! scalars in
ts,dt,dx,dy, & ! scalars in
lfn_in, & ! arrays in
tbound, & ! scalars out
tend & ! arrays out
Expand Down Expand Up @@ -1172,7 +1172,7 @@ subroutine prop_ls( id, & ! for debug
IMTS,IMTE,JMTS,JMTE, & ! memory dims for tend
ids,ide,jds,jde, & ! domain dims - where lfn exists
its,ite,jts,jte, & ! tile dims - where is tend computed
ts+dt,dx,dy, & ! scalars in
ts+dt,dt,dx,dy, & ! scalars in
lfn1, & ! arrays in
tbound2, & ! scalars out
tend & ! arrays out
Expand Down Expand Up @@ -1329,7 +1329,7 @@ subroutine tend_ls( id, &
tims,time,tjms,tjme, & ! memory dims for tend
ids,ide,jds,jde, & ! domain - nodes where lfn defined
ints,inte,jnts,jnte, & ! region - nodes where tend computed
t,dx,dy, & ! scalars in
t,dt,dx,dy, & ! scalars in
lfn, & ! arrays in
tbound, & ! scalars out
tend & ! arrays out
Expand All @@ -1348,7 +1348,7 @@ subroutine tend_ls( id, &
integer,intent(in)::id,lims,lime,ljms,ljme,tims,time,tjms,tjme
integer, intent(in)::ids,ide,jds,jde,ints,inte,jnts,jnte
real,intent(in)::t ! time
real,intent(in)::dx,dy ! mesh step
real,intent(in)::dt,dx,dy ! mesh step
real,dimension(lims:lime,ljms:ljme),intent(inout)::lfn ! level set function
real,intent(out)::tbound ! max allowed time step
real,dimension(tims:time,tjms:tjme),intent(out)::tend ! tendency (rhs of the level set pde)
Expand All @@ -1359,7 +1359,7 @@ subroutine tend_ls( id, &
#endif

!*** local
real:: diffLx,diffLy,diffRx,diffRy, &
real:: te,diffLx,diffLy,diffRx,diffRy, &
diffCx,diffCy,diff2x,diff2y,grad,rr, &
ros_back,ros_wind,ros_slope,advx,advy,scale,nvx,nvy,speed,tanphi
integer::i,j
Expand Down Expand Up @@ -1454,28 +1454,58 @@ subroutine tend_ls( id, &
if(fire_grows_only.gt.0)rr=max(rr,0.)

if(fire_upwind_split.eq.0)then

! get rate of spread
tend(i,j) = -rr*grad ! normal term
if (grad > 0.) then
tbound = max(tbound,rr*(abs(diff2x)/dx+abs(diff2y)/dy)/grad)
endif
te = -rr*grad ! normal term

else

! normal direction backing rate only
te = - ros_back*grad

! advection in wind direction
scale=ros_wind/sqrt(vx(i,j)**2+vy(i,j)**2+eps)
! advection field
advx=vx(i,j)*scale
advy=vy(i,j)*scale
! add advection in slope direction
scale=ros_slope/sqrt(dzfsdx(i,j)**2+dzfsdy(i,j)**2+eps)
! advection field
advx=advx+dzfsdx(i,j)*scale
advy=advy+dzfsdx(i,j)*scale

call crash('prop_ls: fire_upwind_split not done yet')
if (abs(speed)> eps) then
advx=vx(i,j)*ros_wind/speed
advy=vy(i,j)*ros_wind/speed
else
advx=0
advy=0
endif

! advection from slope direction
if(abs(tanphi)>eps) then
advx=advx+dzfsdx(i,j)*ros_slope/tanphi
advy=advy+dzfsdy(i,j)*ros_slope/tanphi
endif

if(fire_upwind_split.eq.1)then

! one-sided upwinding
te = te - max(advx,0.)*diffLx - min(advx,0.)*diffRy &
- max(advy,0.)*diffLy - min(advy,0.)*diffRy


elseif(fire_upwind_split.eq.2)then

! Lax-Friedrichs
call crash('prop_ls: bad fire_upwind_split, Lax-Friedrichs not done yet')

else

call crash('prop_ls: bad fire_upwind_split')

endif
endif

! cfl condition
if (grad > 0.) then
tbound = max(tbound,rr*(abs(diff2x)/dx+abs(diff2y)/dy)/grad)
endif

! add numerical viscosity
tend(i,j)=tend(i,j) + fire_viscosity*abs(rr)*((diffRx-diffLx)+(diffRy-diffLy))
te=te + fire_viscosity*abs(rr)*((diffRx-diffLx)+(diffRy-diffLy))

tend(i,j)=te
#ifdef DEBUG_ANY
rra(i,j)=rr
grada(i,j)=grad
Expand Down
2 changes: 1 addition & 1 deletion wrfv2_fire/test/em_fire/namelist.input
Expand Up @@ -166,7 +166,7 @@
!
! verbosity
fire_print_msg = 1, ! 1 print fire debugging messages
fire_print_file = 1, ! 1 write files for matlab; do not use with parallel runs
fire_print_file = 0, ! 1 write files for matlab; do not use with parallel runs
!
! method selections, do not change from defaults unless you know what you are doing
fire_boundary_guard = -1, ! integer, number of cells to stop when fire close to the domain boundary, -1 turn off
Expand Down

0 comments on commit 4b5d8b2

Please sign in to comment.