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
|
|
|
!!
|
2013-12-10 15:23:28 -02:00
|
|
|
!! Copyright (C) 2008-2013 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
|
|
|
|
!!
|
2013-12-11 10:59:25 -02:00
|
|
|
!! This module provides an interface for temporal integration with
|
|
|
|
!! the stability handling. New integration methods can be added by
|
|
|
|
!! implementing more evolve_* subroutines.
|
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
|
|
|
|
|
2013-12-11 10:59:25 -02:00
|
|
|
! pointer to the temporal integration subroutine
|
|
|
|
!
|
|
|
|
procedure(evolve_euler), pointer, save :: evolve => null()
|
|
|
|
|
2012-07-28 11:47:14 -03:00
|
|
|
! evolution parameters
|
|
|
|
!
|
2013-12-11 22:48:48 -02:00
|
|
|
real , save :: cfl = 0.5d+00
|
2012-07-28 11:47:14 -03:00
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! time variables
|
2012-07-28 11:47:14 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
integer, save :: n = 0
|
2013-12-11 22:48:48 -02:00
|
|
|
real , save :: t = 0.0d+00
|
|
|
|
real , save :: dt = 1.0d+00
|
|
|
|
real , save :: dtn = 1.0d+00
|
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
|
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
public :: initialize_evolution, finalize_evolution
|
2013-12-11 22:48:48 -02:00
|
|
|
public :: advance, new_time_step
|
2012-07-28 11:47:14 -03:00
|
|
|
|
|
|
|
! declare public variables
|
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
public :: cfl, n, t, dt, dtn
|
2012-07-28 11:47:14 -03:00
|
|
|
|
|
|
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
|
!
|
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 *****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-28 11:47:14 -03:00
|
|
|
! subroutine INITIALIZE_EVOLUTION:
|
|
|
|
! -------------------------------
|
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
! Subroutine initializes module EVOLUTION by setting its parameters.
|
2012-07-28 11:47:14 -03:00
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! verbose - a logical flag turning the information printing;
|
|
|
|
! iret - an integer flag for error return value;
|
2012-07-28 11:47:14 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
subroutine initialize_evolution(verbose, iret)
|
2012-07-28 11:47:14 -03:00
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! include external procedures
|
2012-07-28 11:47:14 -03:00
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
use parameters, only : get_parameter_string, get_parameter_real
|
2012-07-28 11:47:14 -03:00
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
2013-12-11 10:59:25 -02:00
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
logical, intent(in) :: verbose
|
|
|
|
integer, intent(inout) :: iret
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
character(len=255) :: integration = "rk2"
|
|
|
|
character(len=255) :: name_int = ""
|
2012-07-28 11:47:14 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
! get the integration method and the value of the CFL coefficient
|
2012-07-28 11:47:14 -03:00
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
call get_parameter_string("time_advance", integration)
|
|
|
|
call get_parameter_real ("cfl" , cfl )
|
2012-07-28 11:47:14 -03:00
|
|
|
|
2013-12-11 10:59:25 -02:00
|
|
|
! select the integration method, check the correctness of the integration
|
|
|
|
! parameters and adjust the CFL coefficient if necessary
|
2012-07-31 15:13:51 -03:00
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
select case(trim(integration))
|
|
|
|
|
|
|
|
case ("euler", "EULER")
|
|
|
|
|
|
|
|
name_int = "1st order Euler method"
|
|
|
|
evolve => evolve_euler
|
|
|
|
|
|
|
|
case ("rk2", "RK2")
|
|
|
|
|
|
|
|
name_int = "2nd order Runge-Kutta method"
|
|
|
|
evolve => evolve_rk2
|
|
|
|
|
|
|
|
case default
|
|
|
|
|
|
|
|
if (verbose) then
|
|
|
|
write (*,"(1x,a)") "The selected time advance method is not " // &
|
|
|
|
"implemented: " // trim(integration)
|
|
|
|
name_int = "2nd order Runge-Kutta method"
|
|
|
|
evolve => evolve_rk2
|
|
|
|
end if
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
! ! calculate the initial time step
|
|
|
|
! !
|
|
|
|
! call new_time_step()
|
|
|
|
|
|
|
|
! print information about the Riemann solver
|
|
|
|
!
|
|
|
|
if (verbose) then
|
|
|
|
|
|
|
|
write (*,"(4x,a,1x,a)" ) "time advance =", trim(name_int)
|
|
|
|
|
|
|
|
end if
|
2012-07-31 15:13:51 -03:00
|
|
|
|
2012-07-28 11:47:14 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine initialize_evolution
|
|
|
|
!
|
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
! subroutine FINALIZE_EVOLUTION:
|
|
|
|
! -----------------------------
|
|
|
|
!
|
|
|
|
! Subroutine releases memory used by the module.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! iret - an integer flag for error return value;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine finalize_evolution(iret)
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer, intent(inout) :: iret
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine finalize_evolution
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-31 15:04:40 -03:00
|
|
|
! subroutine ADVANCE:
|
|
|
|
! ------------------
|
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
! Subroutine advances the solution by one time step using the selected time
|
|
|
|
! integration method.
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine advance()
|
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! include external procedures
|
|
|
|
!
|
|
|
|
use boundaries , only : boundary_variables
|
|
|
|
use mesh , only : update_mesh
|
|
|
|
|
|
|
|
! include external variables
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
use coordinates , only : toplev
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
! find new time step
|
|
|
|
!
|
|
|
|
call new_time_step()
|
|
|
|
|
|
|
|
! advance the solution using the selected method
|
|
|
|
!
|
|
|
|
call evolve()
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! 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
|
|
|
|
|
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
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine advance
|
|
|
|
!
|
|
|
|
!===============================================================================
|
2012-07-28 11:47:14 -03:00
|
|
|
!!
|
|
|
|
!!*** PRIVATE SUBROUTINES ****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
2010-12-01 10:39:18 -02:00
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine EVOLVE_EULER:
|
|
|
|
! -----------------------
|
|
|
|
!
|
|
|
|
! Subroutine advances the solution by one time step using the 1st order
|
|
|
|
! Euler integration method.
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine evolve_euler()
|
|
|
|
|
|
|
|
! include external procedures
|
|
|
|
!
|
|
|
|
use boundaries , only : boundary_variables
|
|
|
|
use schemes , only : update_increment
|
|
|
|
|
|
|
|
! include external variables
|
|
|
|
!
|
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use coordinates , only : adx, ady, adz
|
|
|
|
use coordinates , only : im, jm, km
|
|
|
|
use equations , only : nv
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_data), pointer :: pblock
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension(3) :: dh
|
|
|
|
real, dimension(nv,im,jm,km) :: du
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! update fluxes for the first step of the RK2 integration
|
|
|
|
!
|
|
|
|
call update_fluxes()
|
|
|
|
|
|
|
|
! 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
|
|
|
|
!
|
|
|
|
dh(1) = dt / adx(pblock%meta%level)
|
|
|
|
dh(2) = dt / ady(pblock%meta%level)
|
|
|
|
dh(3) = dt / adz(pblock%meta%level)
|
|
|
|
|
|
|
|
! calculate variable increment for the current block
|
|
|
|
!
|
|
|
|
call update_increment(dh(:), pblock%f(:,:,:,:,:), du(:,:,:,:))
|
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
pblock%u0(1:nv,:,:,:) = pblock%u0(1:nv,:,:,:) + du(1:nv,:,:,:)
|
|
|
|
|
|
|
|
! 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()
|
|
|
|
|
|
|
|
! update primitive variables
|
|
|
|
!
|
|
|
|
call update_variables()
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine evolve_euler
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine EVOLVE_RK2:
|
|
|
|
! ---------------------
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
|
|
|
! Subroutine advances the solution by one time step using the 2nd order
|
|
|
|
! Runge-Kutta time integration method.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
subroutine evolve_rk2()
|
2012-07-31 15:04:40 -03:00
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! include external procedures
|
|
|
|
!
|
|
|
|
use boundaries , only : boundary_variables
|
2012-08-01 16:38:10 -03:00
|
|
|
use schemes , only : update_increment
|
2012-08-01 12:56:52 -03:00
|
|
|
|
|
|
|
! include external variables
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use coordinates , only : adx, ady, adz
|
|
|
|
use coordinates , only : im, jm, km
|
2013-12-10 21:36:35 -02:00
|
|
|
use equations , only : nv
|
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
|
|
|
|
!
|
2013-12-10 21:36:35 -02:00
|
|
|
real, dimension(3) :: dh
|
|
|
|
real, dimension(nv,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-08-01 16:38:10 -03:00
|
|
|
call update_increment(dh(:), pblock%f(:,:,:,:,:), du(:,:,:,:))
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
2013-12-10 21:36:35 -02:00
|
|
|
pblock%u1(1:nv,:,:,:) = pblock%u0(1:nv,:,:,:) + du(1:nv,:,:,:)
|
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-08-01 16:38:10 -03:00
|
|
|
call update_increment(dh(:), pblock%f(:,:,:,:,:), du(:,:,:,:))
|
2012-07-31 15:04:40 -03:00
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
2013-12-10 21:36:35 -02:00
|
|
|
pblock%u0(1:nv,:,:,:) = 0.5d0 * (pblock%u0(1:nv,:,:,:) &
|
|
|
|
+ pblock%u1(1:nv,:,:,:) + du(1:nv,:,:,:))
|
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()
|
|
|
|
|
2013-12-11 10:59:25 -02:00
|
|
|
! update primitive variables
|
|
|
|
!
|
|
|
|
call update_variables()
|
|
|
|
|
2012-07-31 15:04:40 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-11 10:59:25 -02:00
|
|
|
end subroutine evolve_rk2
|
2012-07-31 15:04:40 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
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()
|
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! include external procedures
|
2012-07-31 16:02:59 -03:00
|
|
|
!
|
2012-08-01 17:41:56 -03:00
|
|
|
use boundaries , only : boundary_fluxes
|
2012-08-01 16:38:10 -03:00
|
|
|
use schemes , only : update_flux
|
2012-08-01 12:56:52 -03:00
|
|
|
|
|
|
|
! include external variables
|
|
|
|
!
|
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use coordinates , only : adx, ady, adz
|
2012-07-31 16:02:59 -03:00
|
|
|
|
|
|
|
! 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
|
|
|
|
!
|
2012-08-01 17:41:56 -03:00
|
|
|
call boundary_fluxes()
|
2012-07-31 16:02:59 -03:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
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()
|
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! include external procedures
|
|
|
|
!
|
|
|
|
use equations , only : update_primitive_variables
|
|
|
|
|
|
|
|
! include external variables
|
2012-07-31 16:38:16 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
use blocks , only : block_data, list_data
|
2012-07-31 16:38:16 -03:00
|
|
|
|
|
|
|
! 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-08-01 12:56:52 -03:00
|
|
|
! subroutine NEW_TIME_STEP:
|
|
|
|
! ------------------------
|
2012-07-31 16:23:20 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
! Subroutine estimates the new time step from the maximum speed in the system
|
|
|
|
! and source term constraints.
|
2012-07-31 16:23:20 -03:00
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
subroutine new_time_step()
|
2012-07-31 16:23:20 -03:00
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! include external procedures
|
2012-07-31 16:23:20 -03:00
|
|
|
!
|
2013-12-10 20:56:37 -02:00
|
|
|
use equations , only : maxspeed, cmax, cmax2
|
2012-08-01 12:56:52 -03:00
|
|
|
#ifdef MPI
|
2012-08-03 02:48:38 -03:00
|
|
|
use mpitools , only : reduce_maximum_real, reduce_maximum_integer
|
2012-08-01 12:56:52 -03:00
|
|
|
#endif /* MPI */
|
2012-07-31 16:23:20 -03:00
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! include external variables
|
2012-07-31 16:23:20 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use coordinates , only : adx, ady, adz
|
|
|
|
use coordinates , only : toplev
|
2012-07-31 16:23:20 -03:00
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! local variables are not implicit by default
|
2012-07-31 16:23:20 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
implicit none
|
2012-07-31 16:23:20 -03:00
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! local pointers
|
2010-12-01 10:39:18 -02:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
type(block_data), pointer :: pblock
|
2010-12-01 10:39:18 -02:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer :: iret
|
2011-05-07 09:23:16 -03:00
|
|
|
integer(kind=4) :: lev
|
2012-08-01 12:56:52 -03:00
|
|
|
real :: cm, dx_min
|
2011-05-07 09:23:16 -03:00
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! local parameters
|
2011-05-07 09:23:16 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
real, parameter :: eps = tiny(cmax)
|
2010-12-01 10:39:18 -02:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
! reset the maximum speed, and the highest level
|
2010-12-01 10:39:18 -02:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
cmax = eps
|
2011-05-07 09:23:16 -03:00
|
|
|
lev = 1
|
2008-12-09 14:51:33 -06:00
|
|
|
|
2011-05-07 09:23:16 -03:00
|
|
|
! iterate over all data blocks in order to find the maximum speed among them
|
2012-08-01 12:56:52 -03:00
|
|
|
! and the highest level which is required to obtain the spatial step
|
2011-05-07 09:23:16 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
pblock => list_data
|
|
|
|
do while (associated(pblock))
|
2011-05-07 09:23:16 -03:00
|
|
|
|
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
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
lev = max(lev, pblock%meta%level)
|
2011-05-07 09:23:16 -03:00
|
|
|
|
|
|
|
! obtain the maximum speed for the current block
|
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
cm = maxspeed(pblock%q(:,:,:,:))
|
2008-12-09 14:51:33 -06:00
|
|
|
|
|
|
|
! compare global and local maximum speeds
|
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
cmax = max(cmax, cm)
|
2011-05-07 09:23:16 -03:00
|
|
|
|
|
|
|
! assiociate the pointer with the next block
|
2008-12-09 14:51:33 -06:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
pblock => pblock%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-08-03 02:48:38 -03:00
|
|
|
call reduce_maximum_real (cmax, iret)
|
|
|
|
call reduce_maximum_integer(lev , iret)
|
2010-12-03 18:02:07 -02:00
|
|
|
#endif /* MPI */
|
2011-05-03 00:01:00 -03:00
|
|
|
|
2013-12-10 20:56:37 -02:00
|
|
|
! calculate squared cmax
|
|
|
|
!
|
|
|
|
cmax2 = cmax * cmax
|
|
|
|
|
2012-08-01 12:56:52 -03:00
|
|
|
! find the smallest spatial step
|
|
|
|
!
|
|
|
|
#if NDIMS == 2
|
|
|
|
dx_min = min(adx(lev), ady(lev))
|
|
|
|
#endif /* NDIMS == 2 */
|
|
|
|
#if NDIMS == 3
|
|
|
|
dx_min = min(adx(lev), ady(lev), adz(lev))
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
! calcilate the new time step
|
2011-05-07 09:23:16 -03:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
dtn = dx_min / max(cmax, eps)
|
2011-05-07 09:23:16 -03:00
|
|
|
|
2013-12-11 22:48:48 -02:00
|
|
|
! calculate the new time step
|
|
|
|
!
|
|
|
|
dt = cfl * dtn
|
|
|
|
|
2008-12-08 21:07:10 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2012-08-01 12:56:52 -03:00
|
|
|
end subroutine new_time_step
|
2011-05-19 18:35:36 -03:00
|
|
|
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2012-08-01 16:38:10 -03:00
|
|
|
end module evolution
|