b4step default routines¶
For GeoClaw, see b4step routine in GeoClaw.
Below are the default b4step library routines for Classic and AMRClaw. By default these do nothing. If you wish to specify aux arrays you will need to copy one of these files to your application directory and modify it as needed. Remember to change to Makefile to point to the proper version.
$CLAW/classic/src/1d/b4step1.f90:
1 2 3 4 5 6 7 8 9 10 11 12 13 14 | subroutine b4step1(mbc,mx,meqn,q,xlower,dx,t,dt,maux,aux)
! Called before each call to step1.
! Use to set time-dependent aux arrays or perform other tasks.
!
! This default version does nothing.
implicit none
integer, intent(in) :: mbc,mx,meqn,maux
real(kind=8), intent(in) :: xlower,dx,t,dt
real(kind=8), intent(inout) :: q(meqn,1-mbc:mx+mbc)
real(kind=8), intent(inout) :: aux(maux,1-mbc:mx+mbc)
end subroutine b4step1
|
$CLAW/classic/src/2d/b4step2.f90:
1 2 3 4 5 6 7 8 9 10 11 12 13 14 | subroutine b4step2(mbc,mx,my,meqn,q,xlower,ylower,dx,dy,t,dt,maux,aux)
! Called before each call to step2.
! Use to set time-dependent aux arrays or perform other tasks.
!
! This default version does nothing.
implicit none
integer, intent(in) :: mbc,mx,my,meqn,maux
real(kind=8), intent(in) :: xlower,ylower,dx,dy,t,dt
real(kind=8), intent(inout) :: q(meqn,1-mbc:mx+mbc,1-mbc:my+mbc)
real(kind=8), intent(inout) :: aux(maux,1-mbc:mx+mbc,1-mbc:my+mbc)
end subroutine b4step2
|
$CLAW/classic/src/3d/b4step3.f90:
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 | subroutine b4step3(mbc,mx,my,mz,meqn,q,xlower,ylower,zlower, &
dx,dy,dz,t,dt,maux,aux)
! Called before each call to step3.
! Use to set time-dependent aux arrays or perform other tasks.
!
! This default version does nothing.
implicit none
integer, intent(in) :: mbc,mx,my,mz,meqn,maux
real(kind=8), intent(in) :: xlower,ylower,zlower,dx,dy,dz,t,dt
real(kind=8), intent(inout) :: q(meqn,1-mbc:mx+mbc,1-mbc:my+mbc,1-mbc:mz+mbc)
real(kind=8), intent(inout) :: aux(maux,1-mbc:mx+mbc,1-mbc:my+mbc,1-mbc:mz+mbc)
end subroutine b4step3
|
b4step routine in GeoClaw¶
In GeoClaw, there is a library routine that sets aux(1,i,j) to the cell average of the bathymetry, aux(2,i,j) to the ratio of the true cell area to dx*dy (the capacity function), and aux(3,i,j) to the length ratio of the bottom edge to dx (The latter two quantities vary with latitude when coordinate_system == 2).
$CLAW/geoclaw/src/2d/shallow/b4step2.f90:
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 | ! ============================================
subroutine b4step2(mbc,mx,my,meqn,q,xlower,ylower,dx,dy,t,dt,maux,aux)
! ============================================
!
! # called before each call to step
! # use to set time-dependent aux arrays or perform other tasks.
!
! This particular routine sets negative values of q(1,i,j) to zero,
! as well as the corresponding q(m,i,j) for m=1,meqn.
! This is for problems where q(1,i,j) is a depth.
! This should occur only because of rounding error.
!
! Also calls movetopo if topography might be moving.
use geoclaw_module, only: dry_tolerance
use geoclaw_module, only: g => grav
use topo_module, only: num_dtopo,topotime
use topo_module, only: aux_finalized
use topo_module, only: xlowdtopo,xhidtopo,ylowdtopo,yhidtopo
use amr_module, only: xlowdomain => xlower
use amr_module, only: ylowdomain => ylower
use amr_module, only: xhidomain => xupper
use amr_module, only: yhidomain => yupper
use amr_module, only: xperdom,yperdom,spheredom,NEEDS_TO_BE_SET
use storm_module, only: set_storm_fields
implicit none
! Subroutine arguments
integer, intent(in) :: meqn
integer, intent(inout) :: mbc,mx,my,maux
real(kind=8), intent(inout) :: xlower, ylower, dx, dy, t, dt
real(kind=8), intent(inout) :: q(meqn,1-mbc:mx+mbc,1-mbc:my+mbc)
real(kind=8), intent(inout) :: aux(maux,1-mbc:mx+mbc,1-mbc:my+mbc)
! Local storage
integer :: index,i,j,k,dummy
real(kind=8) :: h,u,v
! Check for NaNs in the solution
call check4nans(meqn,mbc,mx,my,q,t,1)
! check for h < 0 and reset to zero
! check for h < drytolerance
! set hu = hv = 0 in all these cells
forall(i=1-mbc:mx+mbc, j=1-mbc:my+mbc,q(1,i,j) < dry_tolerance)
q(1,i,j) = max(q(1,i,j),0.d0)
q(2:3,i,j) = 0.d0
end forall
if (aux_finalized < 2) then
! topo arrays might have been updated by dtopo more recently than
! aux arrays were set unless at least 1 step taken on all levels
aux(1,:,:) = NEEDS_TO_BE_SET ! new system checks this val before setting
call setaux(mbc,mx,my,xlower,ylower,dx,dy,maux,aux)
endif
! Set wind and pressure aux variables for this grid
call set_storm_fields(maux,mbc,mx,my,xlower,ylower,dx,dy,t,aux)
end subroutine b4step2
|