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 15:13:51 -03:00
|
|
|
public :: initialize_evolution, advance, evolve
|
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 equations , only : update_primitive_variables
|
|
|
|
use mesh , only : update_mesh
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
type(block_data), pointer :: pblock
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
#ifdef RK2
|
|
|
|
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 */
|
|
|
|
|
|
|
|
! update solution using numerical fluxes stored in 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
|
|
|
|
|
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
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
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
|
|
|
!
|
2010-12-01 15:14:07 -02:00
|
|
|
subroutine evolve()
|
2008-12-07 18:57:08 -06:00
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use blocks , only : block_data, list_data
|
|
|
|
use boundaries , only : boundary_variables
|
2011-04-30 15:10:02 -03:00
|
|
|
#ifdef CONSERVATIVE
|
2012-07-27 17:01:21 -03:00
|
|
|
use boundaries , only : boundary_correct_fluxes
|
2011-04-30 15:10:02 -03:00
|
|
|
#endif /* CONSERVATIVE */
|
2011-04-26 12:37:25 -03:00
|
|
|
#ifdef REFINE
|
2012-07-27 17:01:21 -03:00
|
|
|
use coordinates, only : toplev
|
2011-04-26 12:37:25 -03:00
|
|
|
#endif /* REFINE */
|
2012-07-27 22:33:18 -03:00
|
|
|
use equations , only : update_primitive_variables
|
2012-07-27 17:01:21 -03:00
|
|
|
use mesh , only : update_mesh
|
2011-05-07 09:23:16 -03:00
|
|
|
#ifdef FORCE
|
2012-07-27 17:01:21 -03:00
|
|
|
use forcing , only : tbfor
|
|
|
|
use forcing , only : fourier_transform, evolve_forcing
|
|
|
|
use variables , only : idn, imz
|
2011-05-07 09:23:16 -03:00
|
|
|
#endif /* FORCE */
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2009-09-14 18:28:17 -03:00
|
|
|
type(block_data), pointer :: pblock
|
2010-12-01 10:39:18 -02:00
|
|
|
real :: 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
|
|
|
!
|
2010-12-08 18:19:16 -02:00
|
|
|
#ifdef FORCE
|
2011-06-10 15:46:40 -03:00
|
|
|
! perform forcing evolution only when t >= tbfor
|
|
|
|
!
|
|
|
|
if (t .ge. tbfor) then
|
|
|
|
|
2011-03-07 00:08:31 -03:00
|
|
|
! perform the Fourier transform of the velocity field
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock => list_data
|
|
|
|
do while (associated(pblock))
|
2011-03-07 00:08:31 -03:00
|
|
|
|
2011-06-10 15:46:40 -03:00
|
|
|
if (pblock%meta%leaf) &
|
|
|
|
call fourier_transform(pblock%meta%level &
|
|
|
|
, pblock%meta%xmin, pblock%meta%ymin &
|
|
|
|
, pblock%meta%zmin, pblock%u(idn:imz,:,:,:))
|
2011-03-07 00:08:31 -03:00
|
|
|
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock => pblock%next ! assign pointer to the next block
|
2011-03-07 00:08:31 -03:00
|
|
|
|
2011-06-10 15:46:40 -03:00
|
|
|
end do
|
2011-03-07 00:08:31 -03:00
|
|
|
|
2010-12-08 18:19:16 -02:00
|
|
|
! evolve the forcing source terms by the time interval dt
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
call evolve_forcing(t, dt)
|
|
|
|
|
|
|
|
end if ! t >= tbfor
|
2010-12-08 18:19:16 -02:00
|
|
|
#endif /* FORCE */
|
|
|
|
|
2009-09-14 18:28:17 -03:00
|
|
|
! iterate over all data blocks and perform one step of time evolution
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2009-09-14 18:28:17 -03:00
|
|
|
pblock => list_data
|
2008-12-09 14:51:33 -06:00
|
|
|
do while (associated(pblock))
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
! check if this block is a leaf
|
|
|
|
!
|
2011-04-30 12:28:02 -03:00
|
|
|
#ifdef CONSERVATIVE
|
2011-04-30 12:53:17 -03:00
|
|
|
if (pblock%meta%leaf) &
|
|
|
|
#ifdef EULER
|
|
|
|
call flux_euler(pblock)
|
|
|
|
#endif /* EULER */
|
2011-05-01 08:57:14 -03:00
|
|
|
#ifdef RK2
|
|
|
|
call flux_rk2(pblock)
|
|
|
|
#endif /* RK2 */
|
2011-05-01 10:31:27 -03:00
|
|
|
#ifdef RK3
|
|
|
|
call flux_rk3(pblock)
|
|
|
|
#endif /* RK3 */
|
2011-04-30 12:28:02 -03:00
|
|
|
#else /* CONSERVATIVE */
|
2009-09-14 18:28:17 -03:00
|
|
|
if (pblock%meta%leaf) &
|
2008-12-09 21:06:21 -06:00
|
|
|
#ifdef EULER
|
|
|
|
call evolve_euler(pblock)
|
|
|
|
#endif /* EULER */
|
2008-12-08 16:21:59 -06:00
|
|
|
#ifdef RK2
|
2008-12-09 14:51:33 -06:00
|
|
|
call evolve_rk2(pblock)
|
2008-12-08 16:21:59 -06:00
|
|
|
#endif /* RK2 */
|
2010-12-03 15:38:48 -02:00
|
|
|
#ifdef RK3
|
|
|
|
call evolve_rk3(pblock)
|
|
|
|
#endif /* RK3 */
|
2011-04-30 12:28:02 -03:00
|
|
|
#endif /* CONSERVATIVE */
|
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
|
|
|
|
2011-04-30 13:14:07 -03:00
|
|
|
#ifdef CONSERVATIVE
|
2011-04-30 15:10:02 -03:00
|
|
|
! correct the numerical fluxes between neighboring blocks which are at different
|
|
|
|
! levels
|
|
|
|
!
|
2011-05-01 10:44:31 -03:00
|
|
|
call boundary_correct_fluxes()
|
2011-04-30 15:10:02 -03:00
|
|
|
|
2011-04-30 13:14:07 -03:00
|
|
|
! update solution using numerical fluxes stored in data blocks
|
|
|
|
!
|
|
|
|
pblock => list_data
|
|
|
|
do while (associated(pblock))
|
|
|
|
|
|
|
|
! check if this block is a leaf and update its conserved variables using
|
|
|
|
! corrected numerical fluxes
|
|
|
|
!
|
|
|
|
if (pblock%meta%leaf) call update_solution(pblock)
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pblock => pblock%next
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
#endif /* CONSERVATIVE */
|
2008-12-09 20:37:31 -06:00
|
|
|
! update boundaries
|
2008-12-08 12:14:13 -06:00
|
|
|
!
|
2010-12-01 15:14:07 -02:00
|
|
|
call boundary_variables()
|
2008-12-08 12:14:13 -06:00
|
|
|
|
2011-03-10 01:08:00 -03:00
|
|
|
#ifdef REFINE
|
2011-04-26 12:37:25 -03:00
|
|
|
! chec if we need to perform the refinement step
|
|
|
|
!
|
2011-06-06 17:31:51 -03:00
|
|
|
if (toplev .gt. 1) then
|
2011-04-26 12:37:25 -03:00
|
|
|
|
2010-03-14 15:40:24 -03:00
|
|
|
! check refinement and refine
|
|
|
|
!
|
2011-04-26 12:37:25 -03:00
|
|
|
call update_mesh()
|
2010-03-14 15:40:24 -03:00
|
|
|
|
|
|
|
! update boundaries
|
|
|
|
!
|
2011-04-26 12:37:25 -03:00
|
|
|
call boundary_variables()
|
|
|
|
|
2011-06-06 17:31:51 -03:00
|
|
|
end if ! toplev > 1
|
2011-03-10 01:08:00 -03:00
|
|
|
#endif /* REFINE */
|
2011-05-03 00:01:00 -03:00
|
|
|
|
2011-05-07 09:23:16 -03:00
|
|
|
! find new time step
|
2008-12-09 14:51:33 -06:00
|
|
|
!
|
2011-05-07 09:23:16 -03:00
|
|
|
call find_new_timestep()
|
2011-05-03 00:01:00 -03:00
|
|
|
|
2012-07-27 22:33:18 -03:00
|
|
|
! update solution using numerical fluxes stored in 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
|
|
|
|
|
2010-12-01 10:39:18 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine evolve
|
|
|
|
!
|
|
|
|
!===============================================================================
|
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 boundaries , only : boundary_correct_fluxes
|
|
|
|
use coordinates, only : im, jm, km
|
|
|
|
use coordinates, only : adx, ady, adz
|
|
|
|
use scheme , only : update_flux
|
|
|
|
use variables , only : nt, nfl
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_data), pointer :: pblock
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: i, j, k, im1, jm1, km1
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension(3) :: dx
|
|
|
|
real, dimension(nt,im,jm,km) :: du
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
#ifdef CONSERVATIVE
|
|
|
|
! iterate over all data blocks and calculate the first step of
|
|
|
|
! the RK2 integration
|
|
|
|
!
|
|
|
|
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)
|
|
|
|
|
|
|
|
! calculate the flux from U0
|
|
|
|
!
|
|
|
|
do n = 1, NDIMS
|
|
|
|
call update_flux(n, dx(n), pblock%u(:,:,:,:), pblock%f(n,:,:,:,:))
|
|
|
|
end do
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pblock => pblock%next
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! correct the numerical fluxes between neighboring blocks which are at different
|
|
|
|
! levels
|
|
|
|
!
|
|
|
|
call boundary_correct_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
|
|
|
|
!
|
|
|
|
dx(1) = adx(pblock%meta%level)
|
|
|
|
dx(2) = ady(pblock%meta%level)
|
|
|
|
dx(3) = adz(pblock%meta%level)
|
|
|
|
|
|
|
|
! reset the increment array du
|
|
|
|
!
|
|
|
|
du(:,:,:,:) = 0.0d0
|
|
|
|
|
|
|
|
! perform update along the X direction
|
|
|
|
!
|
|
|
|
do i = 2, im
|
|
|
|
im1 = i - 1
|
|
|
|
|
|
|
|
du(:,i,:,:) = du(:,i,:,:) - (pblock%f(1,:,i,:,:) - pblock%f(1,:,im1,:,:)) / dx(1)
|
|
|
|
end do
|
|
|
|
du(:,1,:,:) = du(:,1,:,:) - pblock%f(1,:,1,:,:) / dx(1)
|
|
|
|
|
|
|
|
! perform update along the Y direction
|
|
|
|
!
|
|
|
|
do j = 2, jm
|
|
|
|
jm1 = j - 1
|
|
|
|
|
|
|
|
du(:,:,j,:) = du(:,:,j,:) - (pblock%f(2,:,:,j,:) - pblock%f(2,:,:,jm1,:)) / dx(2)
|
|
|
|
end do
|
|
|
|
du(:,:,1,:) = du(:,:,1,:) - pblock%f(2,:,:,1,:) / dx(2)
|
|
|
|
|
|
|
|
#if NDIMS == 3
|
|
|
|
! perform update along the Z direction
|
|
|
|
!
|
|
|
|
do k = 2, km
|
|
|
|
km1 = k - 1
|
|
|
|
|
|
|
|
du(:,:,:,k) = du(:,:,:,k) - (pblock%f(3,:,:,:,k) - pblock%f(3,:,:,:,km1)) / dx(3)
|
|
|
|
end do
|
|
|
|
du(:,:,:,1) = du(:,:,:,1) - pblock%f(3,:,:,:,1) / dx(3)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
pblock%u1(1:nfl,:,:,:) = pblock%u0(1:nfl,:,:,:) + dt * du(1:nfl,:,:,:)
|
|
|
|
|
|
|
|
! 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()
|
|
|
|
|
|
|
|
! iterate over all data blocks and calculate the second step of
|
|
|
|
! the RK2 integration
|
|
|
|
!
|
|
|
|
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)
|
|
|
|
|
|
|
|
! calculate the flux from U0
|
|
|
|
!
|
|
|
|
do n = 1, NDIMS
|
|
|
|
call update_flux(n, dx(n), pblock%u(:,:,:,:), pblock%f(n,:,:,:,:))
|
|
|
|
end do
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pblock => pblock%next
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! correct the numerical fluxes between neighboring blocks which are at different
|
|
|
|
! levels
|
|
|
|
!
|
|
|
|
call boundary_correct_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
|
|
|
|
!
|
|
|
|
dx(1) = adx(pblock%meta%level)
|
|
|
|
dx(2) = ady(pblock%meta%level)
|
|
|
|
dx(3) = adz(pblock%meta%level)
|
|
|
|
|
|
|
|
! reset the increment array du
|
|
|
|
!
|
|
|
|
du(:,:,:,:) = 0.0d0
|
|
|
|
|
|
|
|
! perform update along the X direction
|
|
|
|
!
|
|
|
|
do i = 2, im
|
|
|
|
im1 = i - 1
|
|
|
|
|
|
|
|
du(:,i,:,:) = du(:,i,:,:) - (pblock%f(1,:,i,:,:) - pblock%f(1,:,im1,:,:)) / dx(1)
|
|
|
|
end do
|
|
|
|
du(:,1,:,:) = du(:,1,:,:) - pblock%f(1,:,1,:,:) / dx(1)
|
|
|
|
|
|
|
|
! perform update along the Y direction
|
|
|
|
!
|
|
|
|
do j = 2, jm
|
|
|
|
jm1 = j - 1
|
|
|
|
|
|
|
|
du(:,:,j,:) = du(:,:,j,:) - (pblock%f(2,:,:,j,:) - pblock%f(2,:,:,jm1,:)) / dx(2)
|
|
|
|
end do
|
|
|
|
du(:,:,1,:) = du(:,:,1,:) - pblock%f(2,:,:,1,:) / dx(2)
|
|
|
|
|
|
|
|
#if NDIMS == 3
|
|
|
|
! perform update along the Z direction
|
|
|
|
!
|
|
|
|
do k = 2, km
|
|
|
|
km1 = k - 1
|
|
|
|
|
|
|
|
du(:,:,:,k) = du(:,:,:,k) - (pblock%f(3,:,:,:,k) - pblock%f(3,:,:,:,km1)) / dx(3)
|
|
|
|
end do
|
|
|
|
du(:,:,:,1) = du(:,:,:,1) - pblock%f(3,:,:,:,1) / dx(3)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
pblock%u0(1:nfl,:,:,:) = 0.5d0 * (pblock%u0(1:nfl,:,:,:) + pblock%u1(1:nfl,:,:,:) + dt * du(1:nfl,:,:,:))
|
|
|
|
|
|
|
|
! 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()
|
|
|
|
#endif /* CONSERVATIVE */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine advance_rk2
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
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
|
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 */
|
|
|
|
use scheme , only : maxspeed, 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
|
|
|
|
!
|
|
|
|
cm = maxspeed(pdata%u(:,:,:,:))
|
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-04-30 12:28:02 -03:00
|
|
|
#ifdef CONSERVATIVE
|
2011-04-30 13:14:07 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! update_solution: subroutine performs an one step update of the conserved
|
|
|
|
! variables for the given data block using the integrated
|
|
|
|
! numerical fluxes stored in the same data block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine update_solution(pblock)
|
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : im, jm, km
|
2012-07-22 22:26:51 -03:00
|
|
|
use coordinates, only : adxi, adyi, adzi
|
2011-05-01 09:44:19 -03:00
|
|
|
#if defined MHD && defined GLM
|
2012-07-27 17:01:21 -03:00
|
|
|
use scheme , only : cmax
|
|
|
|
use variables , only : iph
|
2011-05-01 09:44:19 -03:00
|
|
|
#endif /* MHD & GLM */
|
2011-05-01 09:24:11 -03:00
|
|
|
#ifdef FORCE
|
2012-07-27 17:01:21 -03:00
|
|
|
use forcing , only : tbfor
|
|
|
|
use forcing , only : real_forcing
|
|
|
|
use variables , only : inx, iny, inz
|
|
|
|
use variables , only : idn, imx, imy, imz
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ien
|
2011-05-01 19:21:43 -03:00
|
|
|
#endif /* ADI */
|
2011-05-01 09:24:11 -03:00
|
|
|
#endif /* FORCE */
|
2011-05-07 10:33:49 -03:00
|
|
|
#ifdef SHAPE
|
2012-07-27 19:23:11 -03:00
|
|
|
use problems , only : update_shapes
|
2011-05-07 10:33:49 -03:00
|
|
|
#endif /* SHAPE */
|
2011-04-30 13:14:07 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2011-05-08 00:15:51 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pblock
|
2011-04-30 13:14:07 -03:00
|
|
|
|
2011-05-01 08:57:14 -03:00
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
real :: dxi, dyi, dzi
|
2011-05-01 09:24:11 -03:00
|
|
|
#ifdef FORCE
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension( 3,im,jm,km) :: f
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
real, dimension( im,jm,km) :: ek, dek
|
|
|
|
#endif /* ADI */
|
2011-05-01 09:24:11 -03:00
|
|
|
#endif /* FORCE */
|
2011-05-01 08:57:14 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! prepare dxi, dyi, and dzi
|
|
|
|
!
|
|
|
|
dxi = adxi(pblock%meta%level)
|
|
|
|
dyi = adyi(pblock%meta%level)
|
|
|
|
#if NDIMS == 3
|
|
|
|
dzi = adzi(pblock%meta%level)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
! perform update of conserved variables of the given block using its fluxes
|
|
|
|
!
|
2011-05-19 15:00:20 -03:00
|
|
|
call advance_solution(pblock%u(:,:,:,:), pblock%f(:,:,:,:,:) &
|
|
|
|
, dt, dxi, dyi, dzi)
|
2011-05-01 09:44:19 -03:00
|
|
|
#if defined MHD && defined GLM
|
|
|
|
|
|
|
|
! evolve Psi due to the source term
|
|
|
|
!
|
|
|
|
pblock%u(iph,:,:,:) = decay * pblock%u(iph,:,:,:)
|
|
|
|
#endif /* MHD & GLM */
|
2011-06-10 15:46:40 -03:00
|
|
|
|
2011-05-01 09:24:11 -03:00
|
|
|
#ifdef FORCE
|
2011-06-10 15:46:40 -03:00
|
|
|
! add forcing term only if t >= tbfor
|
|
|
|
!
|
|
|
|
if (t .ge. tbfor) then
|
2011-05-01 09:24:11 -03:00
|
|
|
|
|
|
|
! obtain the forcing terms in real space
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
call real_forcing(pblock%meta%level, pblock%meta%xmin, pblock%meta%ymin &
|
|
|
|
, pblock%meta%zmin, f(:,:,:,:))
|
2011-05-01 09:24:11 -03:00
|
|
|
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
! calculate kinetic energy before adding the forcing term
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
ek(:,:,:) = 0.5d0 * (pblock%u(imx,:,:,:)**2 + pblock%u(imy,:,:,:)**2 &
|
2011-05-01 19:21:43 -03:00
|
|
|
+ pblock%u(imz,:,:,:)**2) / pblock%u(idn,:,:,:)
|
|
|
|
#endif /* ADI */
|
|
|
|
|
2011-05-01 09:24:11 -03:00
|
|
|
! update momenta due to the forcing terms
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(imx,:,:,:) = pblock%u(imx,:,:,:) &
|
2011-05-01 09:24:11 -03:00
|
|
|
+ pblock%u(idn,:,:,:) * f(inx,:,:,:)
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(imy,:,:,:) = pblock%u(imy,:,:,:) &
|
2011-05-01 09:24:11 -03:00
|
|
|
+ pblock%u(idn,:,:,:) * f(iny,:,:,:)
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(imz,:,:,:) = pblock%u(imz,:,:,:) &
|
2011-05-01 09:24:11 -03:00
|
|
|
+ pblock%u(idn,:,:,:) * f(inz,:,:,:)
|
2011-05-01 19:21:43 -03:00
|
|
|
|
|
|
|
#ifdef ADI
|
|
|
|
! calculate kinetic energy after adding the forcing term
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
dek(:,:,:) = 0.5d0 * (pblock%u(imx,:,:,:)**2 + pblock%u(imy,:,:,:)**2 &
|
2011-05-01 19:21:43 -03:00
|
|
|
+ pblock%u(imz,:,:,:)**2) / pblock%u(idn,:,:,:) - ek(:,:,:)
|
|
|
|
|
|
|
|
! update total energy with the injected one
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(ien,:,:,:) = pblock%u(ien,:,:,:) + dek(:,:,:)
|
2011-05-01 19:21:43 -03:00
|
|
|
#endif /* ADI */
|
2011-06-10 15:46:40 -03:00
|
|
|
|
|
|
|
end if ! t >= tbfor
|
2011-05-01 09:24:11 -03:00
|
|
|
#endif /* FORCE */
|
2011-05-01 08:57:14 -03:00
|
|
|
|
2011-05-07 10:33:49 -03:00
|
|
|
#ifdef SHAPE
|
|
|
|
! update solid shapes
|
|
|
|
!
|
2011-05-07 12:14:34 -03:00
|
|
|
call update_shapes(pblock, t)
|
2011-05-07 10:33:49 -03:00
|
|
|
#endif /* SHAPE */
|
|
|
|
|
2011-05-01 08:57:14 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine update_solution
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! advance_solution: subroutine performs an one step update of the conserved
|
|
|
|
! variables array using the numerical fluxes passed as an
|
|
|
|
! argument
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-19 15:00:20 -03:00
|
|
|
subroutine advance_solution(u, f, dh, dxi, dyi, dzi)
|
2011-05-01 08:57:14 -03:00
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use coordinates, only : im, jm, km
|
|
|
|
use variables , only : nqt, nfl
|
|
|
|
use variables , only : inx, iny, inz
|
2011-05-01 09:44:19 -03:00
|
|
|
#ifdef MHD
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ibx, ibz
|
2011-05-01 09:44:19 -03:00
|
|
|
#ifdef GLM
|
2012-07-27 17:01:21 -03:00
|
|
|
use scheme , only : cmax
|
|
|
|
use variables , only : iph
|
2011-05-01 09:44:19 -03:00
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
2011-05-01 08:57:14 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
real, dimension( nqt,im,jm,km), intent(inout) :: u
|
|
|
|
real, dimension(NDIMS,nqt,im,jm,km), intent(in) :: f
|
2011-05-19 15:00:20 -03:00
|
|
|
real :: dh, dxi, dyi, dzi
|
2011-05-01 08:57:14 -03:00
|
|
|
|
2011-04-30 13:14:07 -03:00
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: i, j, k, im1, jm1, km1
|
|
|
|
real :: dhx, dhy, dhz
|
2011-05-01 09:44:19 -03:00
|
|
|
#if defined MHD && defined GLM
|
|
|
|
real :: ch2
|
|
|
|
#endif /* MHD & GLM */
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension(nqt,im,jm,km) :: du
|
2011-04-30 13:14:07 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! prepare dxi, dyi, and dzi
|
|
|
|
!
|
2011-05-19 15:00:20 -03:00
|
|
|
dhx = dh * dxi
|
|
|
|
dhy = dh * dyi
|
2011-04-30 13:14:07 -03:00
|
|
|
#if NDIMS == 3
|
2011-05-19 15:00:20 -03:00
|
|
|
dhz = dh * dzi
|
2011-04-30 13:14:07 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
2011-05-01 09:44:19 -03:00
|
|
|
! reset the increment array du
|
|
|
|
!
|
|
|
|
du(:,:,:,:) = 0.0d0
|
|
|
|
|
2011-04-30 13:14:07 -03:00
|
|
|
! perform update along the X direction
|
|
|
|
!
|
2011-05-18 17:03:24 -03:00
|
|
|
do i = 2, im
|
|
|
|
im1 = i - 1
|
2011-04-30 13:14:07 -03:00
|
|
|
|
2011-05-19 15:00:20 -03:00
|
|
|
du(:,i,:,:) = du(:,i,:,:) - dhx * (f(inx,:,i,:,:) - f(inx,:,im1,:,:))
|
2011-04-30 13:14:07 -03:00
|
|
|
end do
|
2011-05-18 17:03:24 -03:00
|
|
|
du(:,1,:,:) = du(:,1,:,:) - dhx * f(inx,:,1,:,:)
|
2011-04-30 13:14:07 -03:00
|
|
|
|
|
|
|
! perform update along the Y direction
|
|
|
|
!
|
2011-05-18 17:03:24 -03:00
|
|
|
do j = 2, jm
|
|
|
|
jm1 = j - 1
|
2011-04-30 13:14:07 -03:00
|
|
|
|
2011-05-19 15:00:20 -03:00
|
|
|
du(:,:,j,:) = du(:,:,j,:) - dhy * (f(iny,:,:,j,:) - f(iny,:,:,jm1,:))
|
2011-04-30 13:14:07 -03:00
|
|
|
end do
|
2011-05-18 17:03:24 -03:00
|
|
|
du(:,:,1,:) = du(:,:,1,:) - dhy * f(iny,:,:,1,:)
|
2011-04-30 13:14:07 -03:00
|
|
|
|
2011-05-18 17:03:24 -03:00
|
|
|
#if NDIMS == 3
|
2011-04-30 13:14:07 -03:00
|
|
|
! perform update along the Z direction
|
|
|
|
!
|
2011-05-18 17:03:24 -03:00
|
|
|
do k = 2, km
|
|
|
|
km1 = k - 1
|
2011-04-30 13:14:07 -03:00
|
|
|
|
2011-05-19 15:00:20 -03:00
|
|
|
du(:,:,:,k) = du(:,:,:,k) - dhz * (f(inz,:,:,:,k) - f(inz,:,:,:,km1))
|
2011-04-30 13:14:07 -03:00
|
|
|
end do
|
2011-05-18 17:03:24 -03:00
|
|
|
du(:,:,:,1) = du(:,:,:,1) - dhz * f(inz,:,:,:,1)
|
2011-04-30 13:14:07 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
2011-05-01 09:44:19 -03:00
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
u( 1:nfl,:,:,:) = u( 1:nfl,:,:,:) + du( 1:nfl,:,:,:)
|
|
|
|
|
|
|
|
#ifdef MHD
|
|
|
|
! update the solution for the magnetic variables
|
|
|
|
!
|
|
|
|
u(ibx:ibz,:,:,:) = u(ibx:ibz,:,:,:) + du(ibx:ibz,:,:,:)
|
|
|
|
|
|
|
|
#ifdef GLM
|
|
|
|
! calculate c_h^2
|
|
|
|
!
|
|
|
|
ch2 = cmax * cmax
|
|
|
|
|
|
|
|
! update the solution for the scalar potential Psi
|
|
|
|
!
|
|
|
|
u(iph,:,:,:) = u(iph,:,:,:) + ch2 * du(iph,:,:,:)
|
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
|
|
|
|
2011-04-30 13:14:07 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-01 08:57:14 -03:00
|
|
|
end subroutine advance_solution
|
2011-05-19 18:35:36 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! advance_solution_1d: subroutine performs an one step update of the conserved
|
|
|
|
! variables array using the numerical fluxes passed as an
|
|
|
|
! argument along one selected direction only
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine advance_solution_1d(idir, dh, dxi, u, f)
|
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use coordinates, only : im, jm, km
|
|
|
|
use variables , only : nqt, nfl
|
|
|
|
use variables , only : inx, iny, inz
|
2011-05-19 18:35:36 -03:00
|
|
|
#ifdef MHD
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ibx, ibz
|
2011-05-19 18:35:36 -03:00
|
|
|
#ifdef GLM
|
2012-07-27 17:01:21 -03:00
|
|
|
use scheme , only : cmax
|
|
|
|
use variables , only : iph
|
2011-05-19 18:35:36 -03:00
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
integer , intent(in) :: idir
|
|
|
|
real , intent(in) :: dh, dxi
|
|
|
|
real, dimension(nqt,im,jm,km), intent(inout) :: u
|
|
|
|
real, dimension(nqt,im,jm,km), intent(in) :: f
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: i, j, k, im1, jm1, km1
|
|
|
|
real :: dhx, dhy, dhz
|
|
|
|
#if defined MHD && defined GLM
|
|
|
|
real :: ch2
|
|
|
|
#endif /* MHD & GLM */
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension(nqt,im,jm,km) :: du
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! reset the increment array du
|
|
|
|
!
|
|
|
|
du(:,:,:,:) = 0.0d0
|
|
|
|
|
|
|
|
! calculate the conserved variables increment
|
|
|
|
!
|
|
|
|
select case(idir)
|
|
|
|
case(1)
|
|
|
|
|
|
|
|
! prepare dxi
|
|
|
|
!
|
|
|
|
dhx = dh * dxi
|
|
|
|
|
|
|
|
! perform update along the X direction
|
|
|
|
!
|
|
|
|
do i = 2, im
|
|
|
|
im1 = i - 1
|
|
|
|
|
|
|
|
du(:,i,:,:) = du(:,i,:,:) - dhx * (f(:,i,:,:) - f(:,im1,:,:))
|
|
|
|
end do
|
|
|
|
du(:,1,:,:) = du(:,1,:,:) - dhx * f(:,1,:,:)
|
|
|
|
|
|
|
|
case(2)
|
|
|
|
|
|
|
|
! prepare dxi
|
|
|
|
!
|
|
|
|
dhy = dh * dxi
|
|
|
|
|
|
|
|
! perform update along the Y direction
|
|
|
|
!
|
|
|
|
do j = 2, jm
|
|
|
|
jm1 = j - 1
|
|
|
|
|
|
|
|
du(:,:,j,:) = du(:,:,j,:) - dhy * (f(:,:,j,:) - f(:,:,jm1,:))
|
|
|
|
end do
|
|
|
|
du(:,:,1,:) = du(:,:,1,:) - dhy * f(:,:,1,:)
|
|
|
|
|
|
|
|
#if NDIMS == 3
|
|
|
|
case(3)
|
|
|
|
|
|
|
|
! prepare dxi
|
|
|
|
!
|
|
|
|
dhz = dh * dxi
|
|
|
|
|
|
|
|
! perform update along the Z direction
|
|
|
|
!
|
|
|
|
do k = 2, km
|
|
|
|
km1 = k - 1
|
|
|
|
|
|
|
|
du(:,:,:,k) = du(:,:,:,k) - dhz * (f(:,:,:,k) - f(:,:,:,km1))
|
|
|
|
end do
|
|
|
|
du(:,:,:,1) = du(:,:,:,1) - dhz * f(:,:,:,1)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
u( 1:nfl,:,:,:) = u( 1:nfl,:,:,:) + du( 1:nfl,:,:,:)
|
|
|
|
|
|
|
|
#ifdef MHD
|
|
|
|
! update the solution for the magnetic variables
|
|
|
|
!
|
|
|
|
u(ibx:ibz,:,:,:) = u(ibx:ibz,:,:,:) + du(ibx:ibz,:,:,:)
|
|
|
|
|
|
|
|
#ifdef GLM
|
|
|
|
! calculate c_h^2
|
|
|
|
!
|
|
|
|
ch2 = cmax * cmax
|
|
|
|
|
|
|
|
! update the solution for the scalar potential Psi
|
|
|
|
!
|
|
|
|
u(iph,:,:,:) = u(iph,:,:,:) + ch2 * du(iph,:,:,:)
|
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine advance_solution_1d
|
2011-04-30 12:53:17 -03:00
|
|
|
#ifdef EULER
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! flux_euler: subroutine performs the first order integration of the numerical
|
|
|
|
! flux
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine flux_euler(pblock)
|
|
|
|
|
|
|
|
use blocks , only : block_data
|
2012-07-22 22:26:51 -03:00
|
|
|
use coordinates, only : adx, ady, adz
|
2011-04-30 12:53:17 -03:00
|
|
|
use scheme , only : update_flux
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2011-05-08 00:15:51 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pblock
|
2011-04-30 12:53:17 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2011-05-19 18:24:09 -03:00
|
|
|
real :: dx, dy, dz
|
2011-04-30 12:53:17 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! prepare dxi, dyi, and dzi
|
|
|
|
!
|
2011-05-19 18:24:09 -03:00
|
|
|
dx = adx(pblock%meta%level)
|
|
|
|
dy = ady(pblock%meta%level)
|
|
|
|
dz = adz(pblock%meta%level)
|
2011-04-30 12:53:17 -03:00
|
|
|
|
|
|
|
! 1st step of integration
|
|
|
|
!
|
2011-05-19 18:24:09 -03:00
|
|
|
call update_flux(1, dx, pblock%u(:,:,:,:), pblock%f(1,:,:,:,:))
|
|
|
|
call update_flux(2, dy, pblock%u(:,:,:,:), pblock%f(2,:,:,:,:))
|
|
|
|
#if NDIMS == 3
|
|
|
|
call update_flux(3, dz, pblock%u(:,:,:,:), pblock%f(3,:,:,:,:))
|
|
|
|
#endif /* NDIMS == 3 */
|
2011-04-30 12:53:17 -03:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine flux_euler
|
|
|
|
#endif /* EULER */
|
2011-05-01 08:57:14 -03:00
|
|
|
#ifdef RK2
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! flux_rk2: subroutine performs integration of the numerical flux using
|
|
|
|
! the second order Runge-Kutta method
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine flux_rk2(pblock)
|
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : im, jm, km
|
2012-07-22 22:26:51 -03:00
|
|
|
use coordinates, only : adx, ady, adz, adxi, adyi, adzi
|
2012-07-27 17:01:21 -03:00
|
|
|
use scheme , only : update_flux
|
|
|
|
use variables , only : nqt
|
2011-05-01 08:57:14 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2011-05-08 00:15:51 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pblock
|
2011-05-01 08:57:14 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2011-05-19 18:24:09 -03:00
|
|
|
real :: dx, dy, dz, dxi, dyi, dzi
|
2011-05-01 08:57:14 -03:00
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension( nqt,im,jm,km) :: u
|
2011-05-21 22:12:46 -03:00
|
|
|
real, dimension(NDIMS,nqt,im,jm,km) :: f0, f1
|
2011-05-01 08:57:14 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-19 18:24:09 -03:00
|
|
|
! prepare dx, dy, dz, dxi, dyi, and dzi
|
2011-05-01 08:57:14 -03:00
|
|
|
!
|
2011-05-19 18:24:09 -03:00
|
|
|
dx = adx (pblock%meta%level)
|
|
|
|
dy = ady (pblock%meta%level)
|
|
|
|
dz = adz (pblock%meta%level)
|
2011-05-01 08:57:14 -03:00
|
|
|
dxi = adxi(pblock%meta%level)
|
|
|
|
dyi = adyi(pblock%meta%level)
|
|
|
|
dzi = adzi(pblock%meta%level)
|
|
|
|
|
2011-05-19 18:41:51 -03:00
|
|
|
! copy the initial state to the local array u
|
|
|
|
!
|
|
|
|
u(:,:,:,:) = pblock%u(:,:,:,:)
|
|
|
|
|
2011-05-21 22:12:46 -03:00
|
|
|
! calculate fluxes at the moment t
|
2011-05-19 18:41:51 -03:00
|
|
|
!
|
2011-05-21 22:12:46 -03:00
|
|
|
call update_flux(1, dx, u(:,:,:,:), f0(1,:,:,:,:))
|
|
|
|
call update_flux(2, dy, u(:,:,:,:), f0(2,:,:,:,:))
|
|
|
|
#if NDIMS == 3
|
|
|
|
call update_flux(3, dz, u(:,:,:,:), f0(3,:,:,:,:))
|
|
|
|
#endif /* NDIMS == 3 */
|
2011-05-19 18:41:51 -03:00
|
|
|
|
2011-05-21 22:12:46 -03:00
|
|
|
! advance the solution to (t + dt) using the computed fluxes
|
2011-05-19 18:41:51 -03:00
|
|
|
!
|
2011-05-21 22:12:46 -03:00
|
|
|
call advance_solution(u(:,:,:,:), f0(:,:,:,:,:), dt, dxi, dyi, dzi)
|
2011-05-19 18:41:51 -03:00
|
|
|
|
2011-05-21 22:12:46 -03:00
|
|
|
! calculate fluxes at the moment (t + dt/2)
|
2011-05-19 18:41:51 -03:00
|
|
|
!
|
2011-05-21 22:12:46 -03:00
|
|
|
call update_flux(1, dx, u(:,:,:,:), f1(1,:,:,:,:))
|
2011-05-19 18:24:09 -03:00
|
|
|
call update_flux(2, dy, u(:,:,:,:), f1(2,:,:,:,:))
|
|
|
|
#if NDIMS == 3
|
|
|
|
call update_flux(3, dz, u(:,:,:,:), f1(3,:,:,:,:))
|
|
|
|
#endif /* NDIMS == 3 */
|
2011-05-01 08:57:14 -03:00
|
|
|
|
2011-05-21 22:12:46 -03:00
|
|
|
! average the flux at the time )t + dt/2)
|
2011-05-01 08:57:14 -03:00
|
|
|
!
|
|
|
|
pblock%f(:,:,:,:,:) = 0.5d0 * (f0(:,:,:,:,:) + f1(:,:,:,:,:))
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine flux_rk2
|
|
|
|
#endif /* RK2 */
|
2011-05-01 10:31:27 -03:00
|
|
|
#ifdef RK3
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! flux_rk3: subroutine performs integration of the numerical flux using
|
|
|
|
! the third order Runge-Kutta method
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine flux_rk3(pblock)
|
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : im, jm, km
|
2012-07-22 22:26:51 -03:00
|
|
|
use coordinates, only : adx, ady, adz, adxi, adyi, adzi
|
2012-07-27 17:01:21 -03:00
|
|
|
use scheme , only : update_flux
|
|
|
|
use variables , only : nqt
|
2011-05-01 10:31:27 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2011-05-08 00:15:51 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pblock
|
2011-05-01 10:31:27 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
integer :: lev
|
|
|
|
real :: dth, dx, dy, dz, dxi, dyi, dzi
|
2011-05-01 10:31:27 -03:00
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension( nqt,im,jm,km) :: u
|
|
|
|
real, dimension(NDIMS,nqt,im,jm,km) :: f0, f1, f2
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
! obtain the block level
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-19 20:53:38 -03:00
|
|
|
lev = pblock%meta%level
|
2011-05-01 10:31:27 -03:00
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! calculate the half time step
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-19 20:53:38 -03:00
|
|
|
dth = 0.5d0 * dt
|
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! prepare dx, dy, dz, dxi, dyi, and dzi
|
2011-05-19 20:53:38 -03:00
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
dx = adx (lev)
|
|
|
|
dy = ady (lev)
|
|
|
|
dz = adz (lev)
|
|
|
|
dxi = adxi(lev)
|
|
|
|
dyi = adyi(lev)
|
|
|
|
dzi = adzi(lev)
|
2011-05-01 10:31:27 -03:00
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! copy the initial state to the local array u
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
u(:,:,:,:) = pblock%u(:,:,:,:)
|
2011-05-01 10:31:27 -03:00
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! calculate fluxes at time t
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
call update_flux(1, dx, u(:,:,:,:), f0(1,:,:,:,:))
|
|
|
|
call update_flux(2, dy, u(:,:,:,:), f0(2,:,:,:,:))
|
|
|
|
#if NDIMS == 3
|
|
|
|
call update_flux(3, dz, u(:,:,:,:), f0(3,:,:,:,:))
|
|
|
|
#endif /* NDIMS == 3 */
|
2011-05-01 10:31:27 -03:00
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! advance the solution to (t + dt) using the computed fluxes
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
call advance_solution(u(:,:,:,:), f0(:,:,:,:,:), dt, dxi, dyi, dzi)
|
2011-05-01 10:31:27 -03:00
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! calculate fluxes at time (t + dt)
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
call update_flux(1, dx, u(:,:,:,:), f1(1,:,:,:,:))
|
|
|
|
call update_flux(2, dy, u(:,:,:,:), f1(2,:,:,:,:))
|
|
|
|
#if NDIMS == 3
|
|
|
|
call update_flux(3, dz, u(:,:,:,:), f1(3,:,:,:,:))
|
|
|
|
#endif /* NDIMS == 3 */
|
2011-05-01 10:31:27 -03:00
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! average fluxes at the time (t + dt / 2)
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
f2(:,:,:,:,:) = 0.5d0 * (f0(:,:,:,:,:) + f1(:,:,:,:,:))
|
2011-05-19 20:53:38 -03:00
|
|
|
|
|
|
|
! copy the initial state to the local array u
|
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
u(:,:,:,:) = pblock%u(:,:,:,:)
|
2011-05-01 10:31:27 -03:00
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! advance the solution to (t + dt / 2) using the computed fluxes
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
call advance_solution(u(:,:,:,:), f2(:,:,:,:,:), dth, dxi, dyi, dzi)
|
2011-05-01 10:31:27 -03:00
|
|
|
|
2011-05-22 16:49:18 -03:00
|
|
|
! calculate fluxes at time (t + dt / 2)
|
2011-05-01 10:31:27 -03:00
|
|
|
!
|
2011-05-22 16:49:18 -03:00
|
|
|
call update_flux(1, dx, u(:,:,:,:), f2(1,:,:,:,:))
|
|
|
|
call update_flux(2, dy, u(:,:,:,:), f2(2,:,:,:,:))
|
|
|
|
#if NDIMS == 3
|
|
|
|
call update_flux(3, dz, u(:,:,:,:), f2(3,:,:,:,:))
|
|
|
|
#endif /* NDIMS == 3 */
|
2011-05-01 10:31:27 -03:00
|
|
|
|
|
|
|
! calculate the time averaged flux using Gauss formula
|
|
|
|
!
|
|
|
|
pblock%f(:,:,:,:,:) = (f0(:,:,:,:,:) + f1(:,:,:,:,:) &
|
|
|
|
+ 4.0d0 * f2(:,:,:,:,:)) / 6.0d0
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine flux_rk3
|
|
|
|
#endif /* RK3 */
|
2011-04-30 12:28:02 -03:00
|
|
|
#else /* CONSERVATIVE */
|
2008-12-09 21:06:21 -06:00
|
|
|
#ifdef EULER
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! evolve_euler: subroutine evolves the current block using Euler integration
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine evolve_euler(pblock)
|
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : im, jm, km
|
|
|
|
use coordinates, only : adxi, adyi, adzi
|
2010-12-08 21:48:43 -02:00
|
|
|
#ifdef FORCE
|
2012-07-27 17:01:21 -03:00
|
|
|
use forcing , only : tbfor
|
|
|
|
use forcing , only : real_forcing
|
2010-12-08 21:48:43 -02:00
|
|
|
#endif /* FORCE */
|
2008-12-18 12:18:36 -06:00
|
|
|
#ifdef SHAPE
|
2012-07-27 19:23:11 -03:00
|
|
|
use problems , only : update_shapes
|
2008-12-18 12:18:36 -06:00
|
|
|
#endif /* SHAPE */
|
2012-07-27 17:01:21 -03:00
|
|
|
use scheme , only : update, cmax
|
|
|
|
use variables , only : nqt, nfl
|
2010-12-01 10:53:21 -02:00
|
|
|
#ifdef MHD
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ibx, ibz
|
2010-12-01 10:53:21 -02:00
|
|
|
#ifdef GLM
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : iph
|
2010-12-01 10:53:21 -02:00
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
2010-12-08 21:48:43 -02:00
|
|
|
#ifdef FORCE
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : idn, imx, imy, imz
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ien
|
2011-05-01 19:21:43 -03:00
|
|
|
#endif /* ADI */
|
2010-12-08 21:48:43 -02:00
|
|
|
#endif /* FORCE */
|
2008-12-09 21:06:21 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2011-05-08 00:15:51 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pblock
|
2008-12-09 21:06:21 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2010-12-06 16:20:49 -02:00
|
|
|
real :: dxi, dyi, dzi, ch2
|
2008-12-09 21:06:21 -06:00
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
2010-09-19 00:08:20 +02:00
|
|
|
real, dimension(nqt,im,jm,km) :: du
|
2010-12-08 21:48:43 -02:00
|
|
|
#ifdef FORCE
|
|
|
|
real, dimension( 3,im,jm,km) :: f
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
real, dimension( im,jm,km) :: ek, dek
|
|
|
|
#endif /* ADI */
|
2010-12-08 21:48:43 -02:00
|
|
|
#endif /* FORCE */
|
2008-12-09 21:06:21 -06:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2010-12-06 16:20:49 -02:00
|
|
|
! prepare dxi, dyi, and dzi
|
|
|
|
!
|
|
|
|
dxi = adxi(pblock%meta%level)
|
|
|
|
dyi = adyi(pblock%meta%level)
|
|
|
|
dzi = adzi(pblock%meta%level)
|
|
|
|
|
2011-03-06 23:37:51 -03:00
|
|
|
! 1st step of integration
|
|
|
|
!
|
|
|
|
call update(pblock%u(:,:,:,:), du(:,:,:,:), dxi, dyi, dzi)
|
2008-12-18 12:18:36 -06:00
|
|
|
|
2010-12-01 10:53:21 -02:00
|
|
|
! update the solution for the fluid variables
|
2008-12-09 21:06:21 -06:00
|
|
|
!
|
2010-12-01 10:53:21 -02:00
|
|
|
pblock%u(1:nfl,:,:,:) = pblock%u(1:nfl,:,:,:) + dt * du(1:nfl,:,:,:)
|
|
|
|
|
|
|
|
#ifdef MHD
|
|
|
|
! update the solution for the magnetic variables
|
|
|
|
!
|
|
|
|
pblock%u(ibx:ibz,:,:,:) = pblock%u(ibx:ibz,:,:,:) + dt * du(ibx:ibz,:,:,:)
|
|
|
|
|
|
|
|
#ifdef GLM
|
2010-12-01 10:57:40 -02:00
|
|
|
! calculate c_h^2
|
|
|
|
!
|
|
|
|
ch2 = cmax * cmax
|
|
|
|
|
2010-12-01 10:53:21 -02:00
|
|
|
! update the solution for the scalar potential Psi
|
|
|
|
!
|
2010-12-01 10:57:40 -02:00
|
|
|
pblock%u(iph,:,:,:) = pblock%u(iph,:,:,:) + ch2 * dt * du(iph,:,:,:)
|
2010-12-01 11:20:25 -02:00
|
|
|
|
|
|
|
! evolve Psi due to the source term
|
|
|
|
!
|
|
|
|
pblock%u(iph,:,:,:) = decay * pblock%u(iph,:,:,:)
|
2010-12-01 10:53:21 -02:00
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
2011-06-10 15:46:40 -03:00
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
#ifdef FORCE
|
2011-06-10 15:46:40 -03:00
|
|
|
! add the forcing term only if t >= tbfor
|
|
|
|
!
|
|
|
|
if (t .ge. tbfor) then
|
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
! obtain the forcing terms in real space
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
call real_forcing(pblock%meta%level, pblock%meta%xmin, pblock%meta%ymin &
|
2011-04-26 10:36:16 -03:00
|
|
|
, pblock%meta%zmin, f(:,:,:,:))
|
|
|
|
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
! calculate kinetic energy before adding the forcing term
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
ek(:,:,:) = 0.5d0 * (pblock%u(imx,:,:,:)**2 + pblock%u(imy,:,:,:)**2 &
|
2011-05-01 19:21:43 -03:00
|
|
|
+ pblock%u(imz,:,:,:)**2) / pblock%u(idn,:,:,:)
|
|
|
|
#endif /* ADI */
|
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
! update momenta due to the forcing terms
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(imx,:,:,:) = pblock%u(imx,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(1,:,:,:)
|
|
|
|
pblock%u(imy,:,:,:) = pblock%u(imy,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(2,:,:,:)
|
|
|
|
pblock%u(imz,:,:,:) = pblock%u(imz,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(3,:,:,:)
|
2011-04-26 10:36:16 -03:00
|
|
|
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
! calculate kinetic energy after adding the forcing term
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
dek(:,:,:) = 0.5d0 * (pblock%u(imx,:,:,:)**2 + pblock%u(imy,:,:,:)**2 &
|
2011-05-01 19:21:43 -03:00
|
|
|
+ pblock%u(imz,:,:,:)**2) / pblock%u(idn,:,:,:) - ek(:,:,:)
|
|
|
|
|
|
|
|
! update total energy with the injected one
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(ien,:,:,:) = pblock%u(ien,:,:,:) + dek(:,:,:)
|
2011-05-01 19:21:43 -03:00
|
|
|
#endif /* ADI */
|
2011-06-10 15:46:40 -03:00
|
|
|
|
|
|
|
end if ! t >= tbfor
|
2011-04-26 10:36:16 -03:00
|
|
|
#endif /* FORCE */
|
2011-05-07 10:33:49 -03:00
|
|
|
|
|
|
|
#ifdef SHAPE
|
|
|
|
! restrict update in a defined shape
|
2009-10-28 00:12:18 -02:00
|
|
|
!
|
2011-05-07 12:14:34 -03:00
|
|
|
call update_shapes(pblock, t)
|
2011-05-07 10:33:49 -03:00
|
|
|
#endif /* SHAPE */
|
|
|
|
|
2008-12-09 21:06:21 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine evolve_euler
|
|
|
|
#endif /* EULER */
|
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
|
|
|
!
|
2010-12-03 15:38:48 -02:00
|
|
|
! evolve_rk2: subroutine evolves the current block using the 2nd order
|
|
|
|
! Runge-Kutta method
|
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)
|
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : im, jm, km
|
|
|
|
use coordinates, only : adxi, adyi, adzi
|
2010-12-08 21:55:45 -02:00
|
|
|
#ifdef FORCE
|
2012-07-27 17:01:21 -03:00
|
|
|
use forcing , only : tbfor
|
|
|
|
use forcing , only : real_forcing
|
2010-12-08 21:55:45 -02:00
|
|
|
#endif /* FORCE */
|
2008-12-18 12:18:36 -06:00
|
|
|
#ifdef SHAPE
|
2012-07-27 19:23:11 -03:00
|
|
|
use problems , only : update_shapes
|
2008-12-18 12:18:36 -06:00
|
|
|
#endif /* SHAPE */
|
2012-07-27 17:01:21 -03:00
|
|
|
use scheme , only : update, cmax
|
|
|
|
use variables , only : nqt, nfl
|
2010-12-01 10:53:21 -02:00
|
|
|
#ifdef MHD
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ibx, ibz
|
2010-12-01 10:53:21 -02:00
|
|
|
#ifdef GLM
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : iph
|
2010-12-01 10:53:21 -02:00
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
2010-12-08 21:55:45 -02:00
|
|
|
#ifdef FORCE
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : idn, imx, imy, imz
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ien
|
2011-05-01 19:21:43 -03:00
|
|
|
#endif /* ADI */
|
2010-12-08 21:55:45 -02:00
|
|
|
#endif /* FORCE */
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2011-05-08 00:15:51 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pblock
|
2008-12-07 18:57:08 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2010-12-06 16:20:49 -02:00
|
|
|
real :: dxi, dyi, dzi, ch2
|
2008-12-08 16:21:59 -06:00
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
2010-09-19 00:08:20 +02:00
|
|
|
real, dimension(nqt,im,jm,km) :: u1, du
|
2010-12-08 21:55:45 -02:00
|
|
|
#ifdef FORCE
|
|
|
|
real, dimension( 3,im,jm,km) :: f
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
real, dimension( im,jm,km) :: ek, dek
|
|
|
|
#endif /* ADI */
|
2010-12-08 21:55:45 -02:00
|
|
|
#endif /* FORCE */
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
2010-12-06 16:20:49 -02:00
|
|
|
! prepare dxi, dyi, and dzi
|
|
|
|
!
|
|
|
|
dxi = adxi(pblock%meta%level)
|
|
|
|
dyi = adyi(pblock%meta%level)
|
|
|
|
dzi = adzi(pblock%meta%level)
|
|
|
|
|
2010-12-03 15:38:48 -02:00
|
|
|
!! 1st step of integration
|
|
|
|
!!
|
2010-12-06 16:20:49 -02:00
|
|
|
call update(pblock%u(:,:,:,:), du(:,:,:,:), dxi, dyi, dzi)
|
2008-12-18 12:18:36 -06:00
|
|
|
|
2010-12-01 10:53:21 -02:00
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
u1(1:nfl,:,:,:) = pblock%u(1:nfl,:,:,:) + dt * du(1:nfl,:,:,:)
|
|
|
|
|
|
|
|
#ifdef MHD
|
|
|
|
! update the solution for the magnetic variables
|
2008-12-08 16:21:59 -06:00
|
|
|
!
|
2010-12-01 10:53:21 -02:00
|
|
|
u1(ibx:ibz,:,:,:) = pblock%u(ibx:ibz,:,:,:) + dt * du(ibx:ibz,:,:,:)
|
|
|
|
|
|
|
|
#ifdef GLM
|
2011-03-06 23:37:51 -03:00
|
|
|
! calculate c_h^2
|
|
|
|
!
|
|
|
|
ch2 = cmax * cmax
|
|
|
|
|
2010-12-01 10:53:21 -02:00
|
|
|
! update the solution for the scalar potential Psi
|
|
|
|
!
|
2010-12-01 10:57:40 -02:00
|
|
|
u1(iph,:,:,:) = pblock%u(iph,:,:,:) + ch2 * dt * du(iph,:,:,:)
|
2011-04-26 11:50:57 -03:00
|
|
|
|
2010-12-01 10:53:21 -02:00
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
2008-12-08 16:21:59 -06:00
|
|
|
! 2nd step of integration
|
|
|
|
!
|
2010-12-06 16:20:49 -02:00
|
|
|
call update(u1(:,:,:,:), du(:,:,:,:), dxi, dyi, dzi)
|
2008-12-08 16:21:59 -06:00
|
|
|
|
2010-12-01 10:53:21 -02:00
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
pblock%u(1:nfl,:,:,:) = 0.5d0 * (pblock%u(1:nfl,:,:,:) &
|
2010-12-03 15:38:48 -02:00
|
|
|
+ u1(1:nfl,:,:,:) + dt * du(1:nfl,:,:,:))
|
2010-12-01 10:53:21 -02:00
|
|
|
|
|
|
|
#ifdef MHD
|
|
|
|
! update the solution for the magnetic variables
|
|
|
|
!
|
|
|
|
pblock%u(ibx:ibz,:,:,:) = 0.5d0 * (pblock%u(ibx:ibz,:,:,:) &
|
2010-12-03 15:38:48 -02:00
|
|
|
+ u1(ibx:ibz,:,:,:) + dt * du(ibx:ibz,:,:,:))
|
2010-12-01 10:53:21 -02:00
|
|
|
|
|
|
|
#ifdef GLM
|
|
|
|
! update the solution for the scalar potential Psi
|
2008-12-08 16:21:59 -06:00
|
|
|
!
|
2010-12-01 10:53:21 -02:00
|
|
|
pblock%u(iph,:,:,:) = 0.5d0 * (pblock%u(iph,:,:,:) &
|
2010-12-03 15:38:48 -02:00
|
|
|
+ u1(iph,:,:,:) + ch2 * dt * du(iph,:,:,:))
|
2010-12-01 11:20:25 -02:00
|
|
|
|
|
|
|
! evolve Psi due to the source term
|
|
|
|
!
|
|
|
|
pblock%u(iph,:,:,:) = decay * pblock%u(iph,:,:,:)
|
2011-04-26 11:50:57 -03:00
|
|
|
|
2010-12-01 10:53:21 -02:00
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
2011-06-10 15:46:40 -03:00
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
#ifdef FORCE
|
2011-06-10 15:46:40 -03:00
|
|
|
! add the forcing term only if t >= tbfor
|
|
|
|
!
|
|
|
|
if (t .ge. tbfor) then
|
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
! obtain the forcing terms in real space
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
call real_forcing(pblock%meta%level, pblock%meta%xmin, pblock%meta%ymin &
|
2011-04-26 10:36:16 -03:00
|
|
|
, pblock%meta%zmin, f(:,:,:,:))
|
|
|
|
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
! calculate kinetic energy before adding the forcing term
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
ek(:,:,:) = 0.5d0 * (pblock%u(imx,:,:,:)**2 + pblock%u(imy,:,:,:)**2 &
|
2011-05-01 19:21:43 -03:00
|
|
|
+ pblock%u(imz,:,:,:)**2) / pblock%u(idn,:,:,:)
|
|
|
|
#endif /* ADI */
|
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
! update momenta due to the forcing terms
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(imx,:,:,:) = pblock%u(imx,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(1,:,:,:)
|
|
|
|
pblock%u(imy,:,:,:) = pblock%u(imy,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(2,:,:,:)
|
|
|
|
pblock%u(imz,:,:,:) = pblock%u(imz,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(3,:,:,:)
|
2011-04-26 10:36:16 -03:00
|
|
|
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
! calculate kinetic energy after adding the forcing term
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
dek(:,:,:) = 0.5d0 * (pblock%u(imx,:,:,:)**2 + pblock%u(imy,:,:,:)**2 &
|
2011-05-01 19:21:43 -03:00
|
|
|
+ pblock%u(imz,:,:,:)**2) / pblock%u(idn,:,:,:) - ek(:,:,:)
|
|
|
|
|
|
|
|
! update total energy with the injected one
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(ien,:,:,:) = pblock%u(ien,:,:,:) + dek(:,:,:)
|
|
|
|
|
|
|
|
end if ! t >= tbfor
|
2011-05-01 19:21:43 -03:00
|
|
|
#endif /* ADI */
|
2011-04-26 10:36:16 -03:00
|
|
|
#endif /* FORCE */
|
2011-05-07 10:33:49 -03:00
|
|
|
|
|
|
|
#ifdef SHAPE
|
|
|
|
! restrict update in a defined shape
|
2009-10-28 00:12:18 -02:00
|
|
|
!
|
2011-05-07 12:14:34 -03:00
|
|
|
call update_shapes(pblock, t)
|
2011-05-07 10:33:49 -03:00
|
|
|
|
|
|
|
#endif /* SHAPE */
|
|
|
|
|
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 */
|
2010-12-03 15:38:48 -02:00
|
|
|
#ifdef RK3
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! evolve_rk3: subroutine evolves the current block using the 3rd order
|
|
|
|
! Runge-Kutta method
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine evolve_rk3(pblock)
|
|
|
|
|
2012-07-27 17:01:21 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : im, jm, km
|
2010-12-09 17:54:59 -02:00
|
|
|
#ifdef FORCE
|
2012-07-27 17:01:21 -03:00
|
|
|
use forcing , only : tbfor
|
|
|
|
use forcing , only : real_forcing
|
2010-12-09 17:54:59 -02:00
|
|
|
#endif /* FORCE */
|
2012-07-22 22:26:51 -03:00
|
|
|
use coordinates, only : adxi, adyi, adzi
|
2010-12-03 15:38:48 -02:00
|
|
|
#ifdef SHAPE
|
2012-07-27 19:23:11 -03:00
|
|
|
use problems , only : update_shapes
|
2010-12-03 15:38:48 -02:00
|
|
|
#endif /* SHAPE */
|
2012-07-27 17:01:21 -03:00
|
|
|
use scheme , only : update, cmax
|
|
|
|
use variables , only : nqt, nfl
|
2010-12-03 15:38:48 -02:00
|
|
|
#ifdef MHD
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ibx, ibz
|
2010-12-03 15:38:48 -02:00
|
|
|
#ifdef GLM
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : iph
|
2010-12-03 15:38:48 -02:00
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
2010-12-09 17:54:59 -02:00
|
|
|
#ifdef FORCE
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : idn, imx, imy, imz
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
2012-07-27 17:01:21 -03:00
|
|
|
use variables , only : ien
|
2011-05-01 19:21:43 -03:00
|
|
|
#endif /* ADI */
|
2010-12-09 17:54:59 -02:00
|
|
|
#endif /* FORCE */
|
2010-12-03 15:38:48 -02:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2011-05-08 00:15:51 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pblock
|
2010-12-03 15:38:48 -02:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2011-05-01 09:44:19 -03:00
|
|
|
real :: dxi, dyi, dzi
|
2010-12-03 15:38:48 -02:00
|
|
|
#if defined MHD && defined GLM
|
2011-05-27 12:59:09 -03:00
|
|
|
real :: ch2
|
2010-12-03 15:38:48 -02:00
|
|
|
#endif /* MHD & GLM */
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
real, dimension(nqt,im,jm,km) :: u1, du
|
2010-12-09 17:54:59 -02:00
|
|
|
#ifdef FORCE
|
|
|
|
real, dimension( 3,im,jm,km) :: f
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
real, dimension( im,jm,km) :: ek, dek
|
|
|
|
#endif /* ADI */
|
2010-12-09 17:54:59 -02:00
|
|
|
#endif /* FORCE */
|
2010-12-03 15:38:48 -02:00
|
|
|
|
|
|
|
! parameters
|
|
|
|
!
|
|
|
|
real, parameter :: f4 = 1.0d0 / 4.0d0, f3 = 1.0d0 / 3.0d0
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2010-12-06 16:20:49 -02:00
|
|
|
! prepare dxi, dyi, and dzi
|
|
|
|
!
|
|
|
|
dxi = adxi(pblock%meta%level)
|
|
|
|
dyi = adyi(pblock%meta%level)
|
|
|
|
dzi = adzi(pblock%meta%level)
|
|
|
|
|
2010-12-03 15:38:48 -02:00
|
|
|
!! 1st step of integration
|
|
|
|
!!
|
2010-12-06 16:20:49 -02:00
|
|
|
call update(pblock%u(:,:,:,:), du(:,:,:,:), dxi, dyi, dzi)
|
2010-12-03 15:38:48 -02:00
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
u1(1:nfl,:,:,:) = pblock%u(1:nfl,:,:,:) + dt * du(1:nfl,:,:,:)
|
|
|
|
|
|
|
|
#ifdef MHD
|
|
|
|
! update the solution for the magnetic variables
|
|
|
|
!
|
|
|
|
u1(ibx:ibz,:,:,:) = pblock%u(ibx:ibz,:,:,:) + dt * du(ibx:ibz,:,:,:)
|
|
|
|
|
|
|
|
#ifdef GLM
|
2011-03-06 23:37:51 -03:00
|
|
|
! calculate c_h^2
|
|
|
|
!
|
|
|
|
ch2 = cmax * cmax
|
|
|
|
|
2010-12-03 15:38:48 -02:00
|
|
|
! update the solution for the scalar potential Psi
|
|
|
|
!
|
|
|
|
u1(iph,:,:,:) = pblock%u(iph,:,:,:) + ch2 * dt * du(iph,:,:,:)
|
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
|
|
|
|
|
|
|
!! 2nd step of integration
|
|
|
|
!!
|
2010-12-06 16:20:49 -02:00
|
|
|
call update(u1(:,:,:,:), du(:,:,:,:), dxi, dyi, dzi)
|
2010-12-03 15:38:48 -02:00
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
u1(1:nfl,:,:,:) = f4 * (3.0d0 * pblock%u(1:nfl,:,:,:) &
|
|
|
|
+ u1(1:nfl,:,:,:) + dt * du(1:nfl,:,:,:))
|
|
|
|
|
|
|
|
#ifdef MHD
|
|
|
|
! update the solution for the magnetic variables
|
|
|
|
!
|
|
|
|
u1(ibx:ibz,:,:,:) = f4 * (3.0d0 * pblock%u(ibx:ibz,:,:,:) &
|
|
|
|
+ u1(ibx:ibz,:,:,:) + dt * du(ibx:ibz,:,:,:))
|
|
|
|
|
|
|
|
#ifdef GLM
|
|
|
|
! update the solution for the scalar potential Psi
|
|
|
|
!
|
|
|
|
u1(iph,:,:,:) = f4 * (3.0d0 * pblock%u(iph,:,:,:) &
|
|
|
|
+ u1(iph,:,:,:) + ch2 * dt * du(iph,:,:,:))
|
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
|
|
|
|
|
|
|
!! 3rd step of integration
|
|
|
|
!!
|
2010-12-06 16:20:49 -02:00
|
|
|
call update(u1(:,:,:,:), du(:,:,:,:), dxi, dyi, dzi)
|
2010-12-03 15:38:48 -02:00
|
|
|
|
|
|
|
! update the solution for the fluid variables
|
|
|
|
!
|
|
|
|
pblock%u(1:nfl,:,:,:) = f3 * (pblock%u(1:nfl,:,:,:) &
|
|
|
|
+ 2.0d0 * (u1(1:nfl,:,:,:) + dt * du(1:nfl,:,:,:)))
|
|
|
|
|
|
|
|
#ifdef MHD
|
|
|
|
! update the solution for the magnetic variables
|
|
|
|
!
|
|
|
|
pblock%u(ibx:ibz,:,:,:) = f3 * (pblock%u(ibx:ibz,:,:,:) &
|
|
|
|
+ 2.0d0 * (u1(ibx:ibz,:,:,:) + dt * du(ibx:ibz,:,:,:)))
|
|
|
|
|
|
|
|
#ifdef GLM
|
|
|
|
! update the solution for the scalar potential Psi
|
|
|
|
!
|
|
|
|
pblock%u(iph,:,:,:) = f3 * (pblock%u(iph,:,:,:) &
|
|
|
|
+ 2.0d0 * (u1(iph,:,:,:) + ch2 * dt * du(iph,:,:,:)))
|
|
|
|
|
|
|
|
! evolve analytically Psi due to the source term
|
|
|
|
!
|
|
|
|
pblock%u(iph,:,:,:) = decay * pblock%u(iph,:,:,:)
|
|
|
|
#endif /* GLM */
|
|
|
|
#endif /* MHD */
|
2011-05-01 19:21:43 -03:00
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
#ifdef FORCE
|
2011-06-10 15:46:40 -03:00
|
|
|
! add the forcing term only if t >= tbfor
|
|
|
|
!
|
|
|
|
if (t .ge. tbfor) then
|
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
! obtain the forcing terms in real space
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
call real_forcing(pblock%meta%level, pblock%meta%xmin, pblock%meta%ymin &
|
2011-04-26 10:36:16 -03:00
|
|
|
, pblock%meta%zmin, f(:,:,:,:))
|
|
|
|
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
! calculate kinetic energy before adding the forcing term
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
ek(:,:,:) = 0.5d0 * (pblock%u(imx,:,:,:)**2 + pblock%u(imy,:,:,:)**2 &
|
2011-05-01 19:21:43 -03:00
|
|
|
+ pblock%u(imz,:,:,:)**2) / pblock%u(idn,:,:,:)
|
|
|
|
#endif /* ADI */
|
|
|
|
|
2011-04-26 10:36:16 -03:00
|
|
|
! update momenta due to the forcing terms
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(imx,:,:,:) = pblock%u(imx,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(1,:,:,:)
|
|
|
|
pblock%u(imy,:,:,:) = pblock%u(imy,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(2,:,:,:)
|
|
|
|
pblock%u(imz,:,:,:) = pblock%u(imz,:,:,:) &
|
|
|
|
+ pblock%u(idn,:,:,:) * f(3,:,:,:)
|
2011-04-26 10:36:16 -03:00
|
|
|
|
2011-05-01 19:21:43 -03:00
|
|
|
#ifdef ADI
|
|
|
|
! calculate kinetic energy after adding the forcing term
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
dek(:,:,:) = 0.5d0 * (pblock%u(imx,:,:,:)**2 + pblock%u(imy,:,:,:)**2 &
|
2011-05-01 19:21:43 -03:00
|
|
|
+ pblock%u(imz,:,:,:)**2) / pblock%u(idn,:,:,:) - ek(:,:,:)
|
|
|
|
|
|
|
|
! update total energy with the injected one
|
|
|
|
!
|
2011-06-10 15:46:40 -03:00
|
|
|
pblock%u(ien,:,:,:) = pblock%u(ien,:,:,:) + dek(:,:,:)
|
|
|
|
|
|
|
|
end if ! t >= tbfor
|
2011-05-01 19:21:43 -03:00
|
|
|
#endif /* ADI */
|
2011-04-26 10:36:16 -03:00
|
|
|
#endif /* FORCE */
|
2011-05-07 10:33:49 -03:00
|
|
|
|
|
|
|
#ifdef SHAPE
|
|
|
|
! restrict update in a defined shape
|
2010-12-03 15:38:48 -02:00
|
|
|
!
|
2011-05-07 12:14:34 -03:00
|
|
|
call update_shapes(pblock, t)
|
2011-05-07 10:33:49 -03:00
|
|
|
|
|
|
|
#endif /* SHAPE */
|
|
|
|
|
2010-12-03 15:38:48 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine evolve_rk3
|
|
|
|
#endif /* RK3 */
|
2011-04-30 12:28:02 -03:00
|
|
|
#endif /* CONSERVATIVE */
|
2010-02-28 18:35:57 -03:00
|
|
|
!
|
2008-12-08 21:07:10 -06:00
|
|
|
!===============================================================================
|
2008-12-07 18:57:08 -06:00
|
|
|
!
|
|
|
|
end module
|