2008-12-07 18:57:08 -06:00
|
|
|
!!*****************************************************************************
|
|
|
|
!!
|
|
|
|
!! module: evolution - handling the time evolution of the block structure
|
|
|
|
!!
|
|
|
|
!! Copyright (C) 2008 Grzegorz Kowal <kowal@astro.wisc.edu>
|
|
|
|
!!
|
|
|
|
!!*****************************************************************************
|
|
|
|
!!
|
|
|
|
!! This file is part of Godunov-AMR.
|
|
|
|
!!
|
|
|
|
!! Godunov-AMR is free software; you can redistribute it and/or modify
|
|
|
|
!! it under the terms of the GNU General Public License as published by
|
|
|
|
!! the Free Software Foundation; either version 3 of the License, or
|
|
|
|
!! (at your option) any later version.
|
|
|
|
!!
|
|
|
|
!! Godunov-AMR is distributed in the hope that it will be useful,
|
|
|
|
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
!! GNU General Public License for more details.
|
|
|
|
!!
|
|
|
|
!! You should have received a copy of the GNU General Public License
|
|
|
|
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
!!
|
|
|
|
!!*****************************************************************************
|
|
|
|
!!
|
|
|
|
!
|
|
|
|
module evolution
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer, save :: n
|
2008-12-08 16:21:59 -06:00
|
|
|
real , save :: t, dt, dtn
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
contains
|
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
! evolve: subroutine sweeps over all leaf blocks and performs one step time
|
|
|
|
! evolution for each according to the selected integration scheme
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
subroutine evolve
|
|
|
|
|
|
|
|
use blocks, only : block, plist
|
2008-12-09 14:51:33 -06:00
|
|
|
use mesh , only : dx_min
|
|
|
|
use scheme, only : maxspeed
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
type(block), pointer :: pblock
|
2008-12-09 14:51:33 -06:00
|
|
|
real :: cmax, cm
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
! iterate over all blocks and perform one step of time evolution
|
|
|
|
!
|
2008-12-09 14:51:33 -06:00
|
|
|
pblock => plist
|
|
|
|
do while (associated(pblock))
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
! check if this block is a leaf
|
|
|
|
!
|
2008-12-08 16:21:59 -06:00
|
|
|
#ifdef RK2
|
2008-12-09 14:51:33 -06:00
|
|
|
if (pblock%leaf .eq. 'T') &
|
|
|
|
call evolve_rk2(pblock)
|
2008-12-08 16:21:59 -06:00
|
|
|
#endif /* RK2 */
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
2008-12-09 14:51:33 -06:00
|
|
|
pblock => pblock%next
|
2008-12-07 18:57:08 -06:00
|
|
|
|
2008-12-09 14:51:33 -06:00
|
|
|
end do
|
2008-12-07 18:57:08 -06:00
|
|
|
|
2008-12-08 12:14:13 -06:00
|
|
|
! TODO: boundary conditions
|
|
|
|
!
|
|
|
|
|
2008-12-09 14:51:33 -06:00
|
|
|
! reset maximum speed
|
|
|
|
!
|
|
|
|
cmax = 1.0e-8
|
|
|
|
|
|
|
|
! iterate over all blocks in order to find the maximum speed
|
|
|
|
!
|
|
|
|
pblock => plist
|
|
|
|
do while (associated(pblock))
|
|
|
|
|
|
|
|
! check if this block is a leaf
|
|
|
|
!
|
|
|
|
if (pblock%leaf .eq. 'T') &
|
|
|
|
cm = maxspeed(pblock%u)
|
|
|
|
|
|
|
|
! compare global and local maximum speeds
|
|
|
|
!
|
|
|
|
cmax = max(cmax, cm)
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pblock => pblock%next
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! get maximum time step
|
|
|
|
!
|
|
|
|
dtn = dx_min / max(cmax, 1.e-8)
|
|
|
|
|
2008-12-08 21:07:10 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
end subroutine evolve
|
2008-12-08 16:21:59 -06:00
|
|
|
#ifdef RK2
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
! evolve_rk2: subroutine evolves the current block using RK2 integration
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
subroutine evolve_rk2(pblock)
|
|
|
|
|
2008-12-08 16:21:59 -06:00
|
|
|
use blocks, only : block, nvars
|
|
|
|
use config, only : igrids, jgrids, kgrids
|
2008-12-08 20:03:01 -06:00
|
|
|
use mesh , only : adxi, adyi, adzi
|
2008-12-08 19:07:42 -06:00
|
|
|
use scheme, only : update
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2008-12-08 16:21:59 -06:00
|
|
|
type(block), pointer, intent(inout) :: pblock
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2008-12-08 16:21:59 -06:00
|
|
|
integer :: q, i, j, k
|
2008-12-08 20:03:01 -06:00
|
|
|
real :: dxi, dyi, dzi
|
2008-12-08 16:21:59 -06:00
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension(nvars,igrids,jgrids,kgrids) :: u1, du
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 20:03:01 -06:00
|
|
|
! prepare dxi, dyi, and dzi
|
|
|
|
!
|
|
|
|
dxi = adxi(pblock%level)
|
|
|
|
dyi = adyi(pblock%level)
|
|
|
|
dzi = adzi(pblock%level)
|
|
|
|
|
2008-12-08 16:21:59 -06:00
|
|
|
! 1st step of integration
|
|
|
|
!
|
2008-12-08 20:03:01 -06:00
|
|
|
call update(pblock%u, du, dxi, dyi, dzi)
|
2008-12-08 16:21:59 -06:00
|
|
|
|
|
|
|
! update solution
|
|
|
|
!
|
|
|
|
do k = 1, kgrids
|
|
|
|
do j = 1, jgrids
|
|
|
|
do i = 1, igrids
|
|
|
|
do q = 1, nvars
|
|
|
|
u1(q,i,j,k) = pblock%u(q,i,j,k) + dt*du(q,i,j,k)
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
|
|
|
! 2nd step of integration
|
|
|
|
!
|
2008-12-08 20:03:01 -06:00
|
|
|
call update(u1, du, dxi, dyi, dzi)
|
2008-12-08 16:21:59 -06:00
|
|
|
|
|
|
|
! update solution
|
|
|
|
!
|
|
|
|
do k = 1, kgrids
|
|
|
|
do j = 1, jgrids
|
|
|
|
do i = 1, igrids
|
|
|
|
do q = 1, nvars
|
2008-12-08 21:07:10 -06:00
|
|
|
pblock%u(q,i,j,k) = 0.5 * (pblock%u(q,i,j,k) + u1(q,i,j,k) &
|
|
|
|
+ dt*du(q,i,j,k))
|
2008-12-08 16:21:59 -06:00
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
2008-12-07 18:57:08 -06:00
|
|
|
|
2008-12-08 21:07:10 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
end subroutine evolve_rk2
|
2008-12-08 16:21:59 -06:00
|
|
|
#endif /* RK2 */
|
2008-12-07 18:57:08 -06:00
|
|
|
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
end module
|