2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-12-07 18:57:08 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! This file is part of the AMUN source code, a program to perform
|
|
|
|
!! Newtonian or relativistic magnetohydrodynamical simulations on uniform or
|
|
|
|
!! adaptive mesh.
|
2008-12-07 18:57:08 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! Copyright (C) 2008-2012 Grzegorz Kowal <grzegorz@amuncode.org>
|
2008-12-07 18:57:08 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! This program 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.
|
2008-12-07 18:57:08 -06:00
|
|
|
!!
|
2011-04-29 11:21:30 -03:00
|
|
|
!! This program is distributed in the hope that it will be useful,
|
2008-12-07 18:57:08 -06:00
|
|
|
!! 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
|
2012-07-22 12:30:20 -03:00
|
|
|
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
|
2008-12-07 18:57:08 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-12-07 18:57:08 -06:00
|
|
|
!!
|
2012-07-28 11:47:14 -03:00
|
|
|
!! module: EVOLUTION
|
|
|
|
!!
|
|
|
|
!! This module performs the time integration using different methods.
|
2012-07-22 12:30:20 -03:00
|
|
|
!!
|
|
|
|
!!******************************************************************************
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
module evolution
|
|
|
|
|
2012-07-28 11:47:14 -03:00
|
|
|
! module variables are not implicit by default
|
|
|
|
!
|
2008-12-07 18:57:08 -06:00
|
|
|
implicit none
|
|
|
|
|
2012-07-28 11:47:14 -03:00
|
|
|
! evolution parameters
|
|
|
|
!
|
|
|
|
real , save :: cfl = 0.5d0
|
|
|
|
|
|
|
|
#if defined MHD && defined GLM
|
|
|
|
! coefficient controlling the decay of scalar potential Psi
|
|
|
|
!
|
|
|
|
real , save :: alpha_p = 0.5d0
|
|
|
|
real , save :: decay = 1.0d0
|
|
|
|
#endif /* MHD & GLM */
|
|
|
|
|
2008-12-07 18:57:08 -06:00
|
|
|
integer, save :: n
|
2011-05-07 09:32:10 -03:00
|
|
|
real , save :: t, dt, dtn, dxmin
|
2008-12-07 18:57:08 -06:00
|
|
|
|
2012-07-28 11:47:14 -03:00
|
|
|
! by default everything is private
|
|
|
|
!
|
|
|
|
private
|
|
|
|
|
|
|
|
! declare public subroutines
|
|
|
|
!
|
2012-07-31 17:40:34 -03:00
|
|
|
public :: initialize_evolution, advance
|
2012-07-28 11:47:14 -03:00
|
|
|
|
|
|
|
! declare public variables
|
|
|
|
!
|
|
|
|
public :: n, cfl, t, dt, dtn, dxmin
|
|
|
|
|
|
|
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
|
!
|
2008-12-07 18:57:08 -06:00
|
|
|
contains
|
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2012-07-28 11:47:14 -03:00
|
|
|
!!
|
|
|
|
!!*** PUBLIC SUBROUTINES *****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine INITIALIZE_EVOLUTION:
|
|
|
|
! -------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine initializes the EVOLUTION module by setting its parameters.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine initialize_evolution()
|
|
|
|
|
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use parameters, only : get_parameter_real
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! get the value of the stability coefficient
|
|
|
|
!
|
|
|
|
call get_parameter_real("cfl", cfl)
|
|
|
|
|
|
|
|
#if defined MHD && defined GLM
|
|
|
|
! get the value of the stability coefficient
|
|
|
|
!
|
|
|
|
call get_parameter_real("alpha_p", alpha_p)
|
|
|
|
|
|
|
|
! calculate decay factor for magnetic field divergence scalar
|
|
|
|
!
|
|
|
|
decay = exp(- alpha_p * cfl)
|
|
|
|
#endif /* MHD & GLM */
|
|
|
|
|
2012-07-31 15:13:51 -03:00
|
|
|
! calculate the initial time step
|
|
|
|
!
|
|
|
|
call find_new_timestep()
|
|
|
|
|
2012-07-28 11:47:14 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine initialize_evolution
|
|
|
|
!
|
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2012-07-31 15:04:40 -03:00
|
|
|
! subroutine ADVANCE:
|
|
|
|
! ------------------
|
|
|
|
!
|
|
|
|
! Subroutine advances the solution by one time step using the selected
|
|
|
|
! time integration method.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine advance()
|
|
|
|
|
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use boundaries , only : boundary_variables
|
|
|
|
#ifdef REFINE
|
|
|
|
use coordinates, only : toplev
|
|
|
|
#endif /* REFINE */
|
|
|
|
use mesh , only : update_mesh
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
type(block_data), pointer :: pblock
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
#ifdef RK2
|
2012-07-31 16:38:16 -03:00
|
|
|
! advance the solution using the 2nd order Runge-Kutta method
|
2012-07-31 15:04:40 -03:00
|
|
|
call advance_rk2()
|
|
|
|
#endif /* RK2 */
|
|
|
|
|
|
|
|
#ifdef REFINE
|
|
|
|
! chec if we need to perform the refinement step
|
|
|
|
!
|
|
|
|
if (toplev > 1) then
|
|
|
|
|
|
|
|
! check refinement and refine
|
|
|
|
!
|
|
|
|
call update_mesh()
|
|
|
|
|
|
|
|
! update boundaries
|
|
|
|
!
|
|
|
|
call boundary_variables()
|
|
|
|
|
|
|
|
end if ! toplev > 1
|
|
|
|
#endif /* REFINE */
|
|
|
|
|
2012-07-31 16:38:16 -03:00
|
|
|
! update primitive variables
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
2012-07-31 16:38:16 -03:00
|
|
|
call update_variables()
|
2012-07-31 15:04:40 -03:00
|
|
|
|
2012-07-31 15:16:25 -03:00
|
|
|
! find new time step
|
|
|
|
!
|
|
|
|
call find_new_timestep()
|
|
|
|
|
2012-07-31 15:04:40 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine advance
|
|
|
|
!
|
|
|
|
!===============================================================================
|
2012-07-28 11:47:14 -03:00
|
|
|
!!
|
|
|
|
!!*** PRIVATE SUBROUTINES ****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
2010-12-01 10:39:18 -02:00
|
|
|
!
|
2012-07-31 15:04:40 -03:00
|
|
|
! subroutine ADVANCE_RK2:
|
|
|
|
! ----------------------
|
|
|
|
!
|
|
|
|
! Subroutine advances the solution by one time step using the 2nd order
|
|
|
|
! Runge-Kutta time integration method.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine advance_rk2()
|
|
|
|
|
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use boundaries , only : boundary_variables
|
|
|
|
use coordinates, only : im, jm, km
|
|
|
|
use coordinates, only : adx, ady, adz
|
2012-07-31 16:23:20 -03:00
|
|
|
use variables , only : nfl
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_data), pointer :: pblock
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
real, dimension(3) :: dh
|
|
|
|
real, dimension(nfl,im,jm,km) :: du
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-31 16:02:59 -03:00
|
|
|
! update fluxes for the first step of the RK2 integration
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
2012-07-31 16:02:59 -03:00
|
|
|
call update_fluxes()
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! update the solution using numerical fluxes stored in the data blocks
|
|
|
|
!
|
|
|
|
pblock => list_data
|
|
|
|
do while (associated(pblock))
|
|
|
|
|
|
|
|
! obtain dx, dy, and dz for the current block
|
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
dh(1) = dt / adx(pblock%meta%level)
|
|
|
|
dh(2) = dt / ady(pblock%meta%level)
|
|
|
|
dh(3) = dt / adz(pblock%meta%level)
|
2012-07-31 15:04:40 -03:00
|
|
|
|
2012-07-31 16:23:20 -03:00
|
|
|
! calculate variable increment for the current block
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
call update_interval(dh(:), pblock%f(:,:,:,:,:), du(:,:,:,:))
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
pblock%u1(1:nfl,:,:,:) = pblock%u0(1:nfl,:,:,:) + du(1:nfl,:,:,:)
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! update the conservative variable pointer
|
|
|
|
!
|
|
|
|
pblock%u => pblock%u1
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pblock => pblock%next
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! update boundaries
|
|
|
|
!
|
|
|
|
call boundary_variables()
|
|
|
|
|
2012-07-31 16:38:16 -03:00
|
|
|
! update primitive variables
|
|
|
|
!
|
|
|
|
call update_variables()
|
|
|
|
|
2012-07-31 16:02:59 -03:00
|
|
|
! update fluxes for the second step of the RK2 integration
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
2012-07-31 16:02:59 -03:00
|
|
|
call update_fluxes()
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! update the solution using numerical fluxes stored in the data blocks
|
|
|
|
!
|
|
|
|
pblock => list_data
|
|
|
|
do while (associated(pblock))
|
|
|
|
|
|
|
|
! obtain dx, dy, and dz for the current block
|
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
dh(1) = dt / adx(pblock%meta%level)
|
|
|
|
dh(2) = dt / ady(pblock%meta%level)
|
|
|
|
dh(3) = dt / adz(pblock%meta%level)
|
2012-07-31 15:04:40 -03:00
|
|
|
|
2012-07-31 16:23:20 -03:00
|
|
|
! calculate variable increment for the current block
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
call update_interval(dh(:), pblock%f(:,:,:,:,:), du(:,:,:,:))
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
pblock%u0(1:nfl,:,:,:) = 0.5d0 * (pblock%u0(1:nfl,:,:,:) &
|
|
|
|
+ pblock%u1(1:nfl,:,:,:) + du(1:nfl,:,:,:))
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! update the conservative variable pointer
|
|
|
|
!
|
|
|
|
pblock%u => pblock%u0
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pblock => pblock%next
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! update boundaries
|
|
|
|
!
|
|
|
|
call boundary_variables()
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine advance_rk2
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-31 16:02:59 -03:00
|
|
|
! subroutine UPDATE_FLUXES:
|
|
|
|
! ------------------------
|
|
|
|
!
|
|
|
|
! Subroutine iterates over all data blocks and calculates the numerical
|
|
|
|
! fluxes for each block. After the fluxes are updated, they are corrected
|
|
|
|
! for blocks which have neighbours at higher refinement level.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine update_fluxes()
|
|
|
|
|
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use boundaries , only : boundary_correct_fluxes
|
|
|
|
use coordinates, only : adx, ady, adz
|
|
|
|
use scheme , only : update_flux
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
type(block_data), pointer :: pblock
|
2012-07-31 16:02:59 -03:00
|
|
|
|
|
|
|
! local vectors
|
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
real, dimension(3) :: dx
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: n
|
2012-07-31 16:02:59 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! iterate over all data blocks
|
|
|
|
!
|
|
|
|
pblock => list_data
|
|
|
|
do while (associated(pblock))
|
|
|
|
|
|
|
|
! obtain dx, dy, and dz for the current block
|
|
|
|
!
|
|
|
|
dx(1) = adx(pblock%meta%level)
|
|
|
|
dx(2) = ady(pblock%meta%level)
|
|
|
|
dx(3) = adz(pblock%meta%level)
|
|
|
|
|
|
|
|
! update the flux for the current block
|
|
|
|
!
|
|
|
|
do n = 1, NDIMS
|
2012-07-31 16:49:14 -03:00
|
|
|
call update_flux(n, dx(n), pblock%q(:,:,:,:), pblock%f(n,:,:,:,:))
|
2012-07-31 16:02:59 -03:00
|
|
|
end do
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pblock => pblock%next
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! correct the numerical fluxes of the blocks which have neighbours at higher
|
|
|
|
! level
|
|
|
|
!
|
|
|
|
call boundary_correct_fluxes()
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine update_fluxes
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-31 16:38:16 -03:00
|
|
|
! subroutine UPDATE_VARIABLES:
|
|
|
|
! ---------------------------
|
|
|
|
!
|
|
|
|
! Subroutine iterates over all data blocks and converts the conservative
|
|
|
|
! variables to their primitive representation.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine update_variables()
|
|
|
|
|
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use equations , only : update_primitive_variables
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_data), pointer :: pblock
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! iterate over all data blocks
|
|
|
|
!
|
|
|
|
pblock => list_data
|
|
|
|
do while (associated(pblock))
|
|
|
|
|
|
|
|
! convert conserved variables to primitive ones for the current block
|
|
|
|
!
|
|
|
|
call update_primitive_variables(pblock%u, pblock%q)
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pblock => pblock%next
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine update_variables
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-31 16:23:20 -03:00
|
|
|
! subroutine UPDATE_INTERVAL:
|
|
|
|
! --------------------------
|
|
|
|
!
|
|
|
|
! Subroutine calculate the conservative variable interval from fluxes.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine update_interval(dh, f, du)
|
|
|
|
|
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use coordinates, only : im, jm, km
|
|
|
|
use variables , only : nfl
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
real, dimension(3) , intent(in) :: dh
|
|
|
|
real, dimension(NDIMS,nfl,im,jm,km), intent(in) :: f
|
|
|
|
real, dimension( nfl,im,jm,km), intent(inout) :: du
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: i, j, k
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! reset the increment array du
|
|
|
|
!
|
|
|
|
du(:,:,:,:) = 0.0d0
|
|
|
|
|
|
|
|
! perform update along the X direction
|
|
|
|
!
|
|
|
|
do i = 2, im
|
|
|
|
du(:,i,:,:) = du(:,i,:,:) - dh(1) * (f(1,:,i,:,:) - f(1,:,i-1,:,:))
|
|
|
|
end do
|
|
|
|
du(:,1,:,:) = du(:,1,:,:) - dh(1) * f(1,:,1,:,:)
|
|
|
|
|
|
|
|
! perform update along the Y direction
|
|
|
|
!
|
|
|
|
do j = 2, jm
|
|
|
|
du(:,:,j,:) = du(:,:,j,:) - dh(2) * (f(2,:,:,j,:) - f(2,:,:,j-1,:))
|
|
|
|
end do
|
|
|
|
du(:,:,1,:) = du(:,:,1,:) - dh(2) * f(2,:,:,1,:)
|
|
|
|
|
|
|
|
#if NDIMS == 3
|
|
|
|
! perform update along the Z direction
|
|
|
|
!
|
|
|
|
do k = 2, km
|
|
|
|
du(:,:,:,k) = du(:,:,:,k) - dh(3) * (f(3,:,:,:,k) - f(3,:,:,:,k-1))
|
|
|
|
end do
|
|
|
|
du(:,:,:,1) = du(:,:,:,1) - dh(3) * f(3,:,:,:,1)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine update_interval
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
! find_new_timestep: subroutine updates the maximum speed among the leafs and
|
|
|
|
! calculates new time step
|
2010-12-01 10:39:18 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
subroutine find_new_timestep()
|
2010-12-01 10:39:18 -02:00
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use blocks , only : block_meta, block_data, list_meta, list_data
|
|
|
|
use coordinates, only : toplev
|
|
|
|
use coordinates, only : adx, ady, adz
|
2012-08-01 12:16:38 -03:00
|
|
|
use equations , only : maxspeed
|
2010-12-03 18:02:07 -02:00
|
|
|
#ifdef MPI
|
2012-07-27 17:01:21 -03:00
|
|
|
use mpitools , only : reduce_maximum_real
|
2010-12-03 18:02:07 -02:00
|
|
|
#endif /* MPI */
|
2012-08-01 12:16:38 -03:00
|
|
|
use variables , only : cmax
|
2010-12-01 10:39:18 -02:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer :: iret
|
2011-05-07 09:32:10 -03:00
|
|
|
real :: cm
|
2011-05-07 09:23:16 -03:00
|
|
|
integer(kind=4) :: lev
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_meta), pointer :: pmeta
|
|
|
|
type(block_data), pointer :: pdata
|
2010-12-01 10:39:18 -02:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
! reset the maximum speed, and highest level
|
2010-12-01 10:39:18 -02:00
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
cmax = 1.0d-16
|
|
|
|
lev = 1
|
2008-12-09 14:51:33 -06:00
|
|
|
|
2011-06-06 17:31:51 -03:00
|
|
|
! if toplev > 1, find the highest level
|
2011-05-07 09:23:16 -03:00
|
|
|
!
|
2011-06-06 17:31:51 -03:00
|
|
|
if (toplev .gt. 1) then
|
2011-05-07 09:23:16 -03:00
|
|
|
|
|
|
|
! iterate over all meta blocks and find the highest level with leafs
|
|
|
|
!
|
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
|
|
|
|
|
|
|
! check if the metablock is a leaf, if so obtaind the highest level
|
|
|
|
!
|
|
|
|
if (pmeta%leaf) lev = max(lev, pmeta%level)
|
|
|
|
|
|
|
|
! associate the pointer with the next block
|
|
|
|
!
|
|
|
|
pmeta => pmeta%next
|
|
|
|
|
|
|
|
end do ! meta blocks
|
|
|
|
|
2011-06-06 17:31:51 -03:00
|
|
|
end if ! toplev > 1
|
2011-05-17 20:26:33 -03:00
|
|
|
|
2011-05-07 09:23:16 -03:00
|
|
|
! find the smallest spacial step
|
|
|
|
!
|
|
|
|
#if NDIMS == 2
|
2011-05-17 20:26:33 -03:00
|
|
|
dxmin = min(adx(lev), ady(lev))
|
2011-05-07 09:23:16 -03:00
|
|
|
#endif /* NDIMS == 2 */
|
|
|
|
#if NDIMS == 3
|
2011-05-17 20:26:33 -03:00
|
|
|
dxmin = min(adx(lev), ady(lev), adz(lev))
|
2011-05-07 09:23:16 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
! iterate over all data blocks in order to find the maximum speed among them
|
|
|
|
!
|
|
|
|
pdata => list_data
|
|
|
|
do while (associated(pdata))
|
2008-12-09 14:51:33 -06:00
|
|
|
|
|
|
|
! check if this block is a leaf
|
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
if (pdata%meta%leaf) then
|
|
|
|
|
2011-06-06 17:31:51 -03:00
|
|
|
! find the maximum level occupied by blocks (can be smaller than toplev)
|
2011-05-07 09:23:16 -03:00
|
|
|
!
|
|
|
|
|
|
|
|
|
|
|
|
! obtain the maximum speed for the current block
|
|
|
|
!
|
2012-07-31 16:55:00 -03:00
|
|
|
cm = maxspeed(pdata%q(:,:,:,:))
|
2008-12-09 14:51:33 -06:00
|
|
|
|
|
|
|
! compare global and local maximum speeds
|
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
cmax = max(cmax, cm)
|
2008-12-09 14:51:33 -06:00
|
|
|
|
2011-05-07 09:23:16 -03:00
|
|
|
end if ! leaf
|
|
|
|
|
|
|
|
! assiociate the pointer with the next block
|
2008-12-09 14:51:33 -06:00
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
pdata => pdata%next
|
2008-12-09 14:51:33 -06:00
|
|
|
|
|
|
|
end do
|
2010-12-03 18:02:07 -02:00
|
|
|
|
|
|
|
#ifdef MPI
|
2012-07-22 19:01:27 -03:00
|
|
|
! find maximum speed in the system from all processors
|
2010-12-03 18:02:07 -02:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
call reduce_maximum_real(cmax, iret)
|
2010-12-03 18:02:07 -02:00
|
|
|
#endif /* MPI */
|
2011-05-03 00:01:00 -03:00
|
|
|
|
2011-05-07 09:23:16 -03:00
|
|
|
! calcilate new time step
|
|
|
|
!
|
|
|
|
dtn = dxmin / max(cmax, 1.0d-16)
|
|
|
|
|
2008-12-08 21:07:10 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
end subroutine find_new_timestep
|
2011-05-19 18:35:36 -03:00
|
|
|
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
end module
|