2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-11 16:12:26 -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-11-11 16:12:26 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! Copyright (C) 2008-2012 Grzegorz Kowal <grzegorz@amuncode.org>
|
2008-11-11 16:12:26 -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-11-11 16:12:26 -06:00
|
|
|
!!
|
2011-04-29 11:21:30 -03:00
|
|
|
!! This program is distributed in the hope that it will be useful,
|
2008-11-11 16:12:26 -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-11-11 16:12:26 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
2012-07-27 19:23:11 -03:00
|
|
|
!! module: PROBLEMS
|
|
|
|
!!
|
|
|
|
!! This module handles the initialization of various test and research
|
|
|
|
!! problems.
|
2012-07-22 12:30:20 -03:00
|
|
|
!!
|
2012-07-27 21:13:43 -03:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!!******************************************************************************
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2012-07-27 19:23:11 -03:00
|
|
|
module problems
|
2008-11-11 16:12:26 -06:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! module variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
2011-05-05 16:51:35 -03:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! module variable to store the problem name
|
|
|
|
!
|
|
|
|
character(len=32), save :: problem = "blast"
|
2011-05-05 16:51:35 -03:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! by default everything is private
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
private
|
|
|
|
|
|
|
|
! declare public subroutines
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
2012-07-27 21:37:57 -03:00
|
|
|
public :: initialize_problems, setup_problem
|
2008-12-22 15:34:02 -06:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
contains
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
2012-07-27 21:13:43 -03:00
|
|
|
!!
|
|
|
|
!!*** PUBLIC SUBROUTINES *****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine INITIALIZE_PROBLEMS:
|
|
|
|
! ------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine prepares module PROBLEMS.
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 10:20:15 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
subroutine initialize_problems()
|
2008-12-18 10:20:15 -06:00
|
|
|
|
2012-07-27 21:28:59 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
2012-07-27 21:42:10 -03:00
|
|
|
use parameters , only : get_parameter_string
|
2012-07-27 21:28:59 -03:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! local variables are not implicit by default
|
2008-12-18 10:20:15 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
implicit none
|
2008-12-18 10:20:15 -06:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
! get the problem name
|
|
|
|
!
|
|
|
|
call get_parameter_string("problem", problem)
|
|
|
|
|
2008-12-18 10:20:15 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
end subroutine initialize_problems
|
2008-12-18 12:18:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
! subroutine SETUP_PROBLEM:
|
|
|
|
! ------------------------
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
! Subroutine sets the initial conditions for selected problem.
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
! Arguments:
|
2011-05-05 16:51:35 -03:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
! pdata - pointer to the datablock structure of the currently initialized
|
|
|
|
! block;
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
2010-12-21 20:25:17 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine setup_problem(pdata)
|
2008-12-22 15:34:02 -06:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! include external procedures and variables
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
2012-07-27 21:42:10 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use error , only : print_error
|
2010-12-21 20:25:17 -02:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! local variables are not implicit by default
|
2010-12-21 20:25:17 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
implicit none
|
2010-12-21 20:25:17 -02:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! input arguments
|
2010-12-21 20:25:17 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pdata
|
2008-12-22 15:34:02 -06:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
! select the setup subroutine depending on the problem name
|
2011-02-16 11:28:50 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
select case(problem)
|
2011-02-16 11:28:50 -02:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! general test problems
|
2010-12-21 20:25:17 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
case("blast")
|
|
|
|
call setup_problem_blast(pdata)
|
2010-12-21 20:25:17 -02:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
case default
|
|
|
|
call print_error("problems::init_problem()" &
|
|
|
|
, "Setup subroutime is not implemented for this problem!")
|
|
|
|
end select
|
2010-12-21 20:25:17 -02:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2010-12-21 20:25:17 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
end subroutine setup_problem
|
2010-12-21 20:25:17 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
!===============================================================================
|
|
|
|
!!
|
|
|
|
!!*** PRIVATE SUBROUTINES ****************************************************
|
|
|
|
!!
|
2011-03-18 14:00:27 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
! subroutine SETUP_PROBLEM_BLAST:
|
|
|
|
! ------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine sets the initial conditions for the spherical blast problem.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pdata - pointer to the datablock structure of the currently initialized
|
|
|
|
! block;
|
|
|
|
!
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2008-12-18 10:20:15 -06:00
|
|
|
!===============================================================================
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
subroutine setup_problem_blast(pdata)
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : im, jm, km
|
2012-07-28 14:11:07 -03:00
|
|
|
use coordinates, only : ax, ay, az, adx, ady, adz
|
2012-07-27 21:42:10 -03:00
|
|
|
use equations , only : prim2cons
|
2012-07-27 21:13:43 -03:00
|
|
|
use equations , only : gamma
|
2012-07-27 21:28:59 -03:00
|
|
|
use parameters , only : get_parameter_real
|
2012-07-27 22:13:52 -03:00
|
|
|
use variables , only : nt
|
|
|
|
use variables , only : idn, ivx, ivy, ivz
|
2010-03-30 23:23:49 -03:00
|
|
|
#ifdef ADI
|
2012-07-27 21:13:43 -03:00
|
|
|
use variables , only : ipr
|
2010-03-30 23:23:49 -03:00
|
|
|
#endif /* ADI */
|
2010-12-01 15:57:58 -02:00
|
|
|
#ifdef MHD
|
2012-07-27 21:13:43 -03:00
|
|
|
use variables , only : ibx, iby, ibz
|
|
|
|
#ifdef GLM
|
|
|
|
use variables , only : iph
|
|
|
|
#endif /* GLM */
|
2010-12-01 15:57:58 -02:00
|
|
|
#endif /* MHD */
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
2008-11-29 22:19:02 -06:00
|
|
|
! input arguments
|
|
|
|
!
|
2011-05-12 09:05:51 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pdata
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! default parameter values
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
real , save :: dens = 1.0d0, ratio = 1.0e2, radius = 0.1d0
|
|
|
|
#ifdef ADI
|
|
|
|
real , save :: csnd = 0.40824829046386301635d0
|
|
|
|
#endif /* ADI */
|
|
|
|
logical, save :: first = .true.
|
2012-07-28 14:11:07 -03:00
|
|
|
real , save :: dn_amb, dn_ovr
|
|
|
|
#ifdef ADI
|
|
|
|
real , save :: pr_amb, pr_ovr
|
|
|
|
#endif /* ADI */
|
|
|
|
real , save :: rad
|
2008-11-29 22:19:02 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-28 14:11:07 -03:00
|
|
|
integer :: i, j, k
|
|
|
|
real :: xl, yl, zl, xu, yu, zu, rl, ru
|
|
|
|
real :: xb, yb, xt, yt
|
|
|
|
real :: dx, dy, dz, dxh, dyh, dzh, daxy
|
|
|
|
real :: fc_amb, fc_ovr
|
2009-09-29 15:56:25 -03:00
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
2012-07-27 22:13:52 -03:00
|
|
|
real, dimension(nt,im) :: q, u
|
|
|
|
real, dimension(im) :: x
|
|
|
|
real, dimension(jm) :: y
|
|
|
|
real, dimension(km) :: z
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2008-12-18 10:20:15 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2010-12-01 15:57:58 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
! prepare problem constants during the first subroutine call
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
if (first) then
|
2010-12-01 15:57:58 -02:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! get problem parameters
|
|
|
|
!
|
|
|
|
call get_parameter_real("dens" , dens )
|
|
|
|
call get_parameter_real("ratio" , ratio )
|
|
|
|
call get_parameter_real("radius", radius)
|
|
|
|
#ifdef ADI
|
|
|
|
call get_parameter_real("csnd" , csnd )
|
2012-07-28 14:11:07 -03:00
|
|
|
#endif /* ADI */
|
|
|
|
|
|
|
|
! calculate the overdense and ambient region densities
|
|
|
|
!
|
|
|
|
dn_amb = dens
|
|
|
|
#ifdef ISO
|
|
|
|
dn_ovr = dn_amb * ratio
|
|
|
|
#endif /* ISO */
|
|
|
|
#ifdef ADI
|
|
|
|
dn_ovr = dn_amb
|
|
|
|
#endif /* ADI */
|
2010-12-01 15:57:58 -02:00
|
|
|
|
2012-07-28 14:11:07 -03:00
|
|
|
#ifdef ADI
|
|
|
|
! calculate parallel and perpendicular pressures from sound speeds
|
2012-07-27 21:13:43 -03:00
|
|
|
!
|
2012-07-28 14:11:07 -03:00
|
|
|
pr_amb = dens * csnd * csnd / gamma
|
|
|
|
pr_ovr = pr_amb * ratio
|
2010-03-30 23:23:49 -03:00
|
|
|
#endif /* ADI */
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2012-07-28 14:11:07 -03:00
|
|
|
! calculate the square of radius
|
|
|
|
!
|
|
|
|
rad = radius * radius
|
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! reset the first execution flag
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
first = .false.
|
2010-12-01 15:57:58 -02:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
end if
|
2010-02-11 23:30:46 -02:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! obtain block coordinates
|
2010-12-01 15:57:58 -02:00
|
|
|
!
|
2012-07-27 22:13:52 -03:00
|
|
|
x(1:im) = pdata%meta%xmin + ax(pdata%meta%level,1:im)
|
|
|
|
y(1:jm) = pdata%meta%ymin + ay(pdata%meta%level,1:jm)
|
2010-12-02 10:02:52 -02:00
|
|
|
#if NDIMS == 3
|
2012-07-27 22:13:52 -03:00
|
|
|
z(1:km) = pdata%meta%zmin + az(pdata%meta%level,1:km)
|
2012-07-27 21:13:43 -03:00
|
|
|
#else /* NDIMS == 3 */
|
2012-07-27 22:13:52 -03:00
|
|
|
z(1:km) = 0.0d0
|
2010-12-02 10:02:52 -02:00
|
|
|
#endif /* NDIMS == 3 */
|
2008-12-16 13:40:34 -06:00
|
|
|
|
2012-07-28 14:11:07 -03:00
|
|
|
! calculate mesh intervals and areas
|
2010-12-01 15:57:58 -02:00
|
|
|
!
|
2012-07-28 14:11:07 -03:00
|
|
|
dx = adx(pdata%meta%level)
|
|
|
|
dy = ady(pdata%meta%level)
|
|
|
|
dz = adz(pdata%meta%level)
|
|
|
|
dxh = 0.5d0 * dx
|
|
|
|
dyh = 0.5d0 * dy
|
|
|
|
#if NDIMS == 3
|
|
|
|
dzh = 0.5d0 * dz
|
|
|
|
#else /* NDIMS == 3 */
|
|
|
|
dzh = 1.0d0
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
daxy = dx * dy
|
|
|
|
|
|
|
|
! set the uniform primitive variables
|
|
|
|
!
|
|
|
|
q(ivx,:) = 0.0d0
|
|
|
|
q(ivy,:) = 0.0d0
|
|
|
|
q(ivz,:) = 0.0d0
|
|
|
|
|
2010-12-01 15:57:58 -02:00
|
|
|
#ifdef MHD
|
2012-07-28 14:11:07 -03:00
|
|
|
! set the uniform magnetic field
|
|
|
|
!
|
|
|
|
q(ibx,:) = bext(inx)
|
|
|
|
q(iby,:) = bext(iny)
|
|
|
|
q(ibz,:) = bext(inz)
|
2012-07-27 21:13:43 -03:00
|
|
|
#endif /* MHD */
|
2010-12-01 15:57:58 -02:00
|
|
|
|
2012-07-28 14:11:07 -03:00
|
|
|
! iterate over all positions in the YZ plane
|
2010-12-01 15:57:58 -02:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
do k = 1, km
|
2012-07-28 14:11:07 -03:00
|
|
|
|
|
|
|
#ifdef R3D
|
|
|
|
! calculate the corner Z coordinates
|
|
|
|
!
|
|
|
|
zl = abs(z(k)) - dzh
|
|
|
|
zu = abs(z(k)) + dzh
|
|
|
|
#endif /* R3D */
|
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
do j = 1, jm
|
2012-07-27 22:13:52 -03:00
|
|
|
|
2012-07-28 14:11:07 -03:00
|
|
|
! calculate the corner Y coordinates
|
|
|
|
!
|
|
|
|
yl = abs(y(j)) - dyh
|
|
|
|
yu = abs(y(j)) + dyh
|
|
|
|
|
|
|
|
! sweep along the X coordinate
|
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
do i = 1, im
|
2010-12-01 15:57:58 -02:00
|
|
|
|
2012-07-28 14:11:07 -03:00
|
|
|
! calculate the corner X coordinates
|
|
|
|
!
|
|
|
|
xl = abs(x(i)) - dxh
|
|
|
|
xu = abs(x(i)) + dxh
|
|
|
|
|
|
|
|
! calculate the minimum and maximum corner distances from the origin
|
2012-07-27 22:13:52 -03:00
|
|
|
!
|
2012-07-28 14:11:07 -03:00
|
|
|
#ifdef R3D
|
|
|
|
rl = xl * xl + yl * yl + zl * zl
|
|
|
|
ru = xu * xu + yu * yu + zu * zu
|
|
|
|
#else /* R3D */
|
|
|
|
rl = xl * xl + yl * yl
|
|
|
|
ru = xu * xu + yu * yu
|
|
|
|
#endif /* R3D */
|
2010-12-01 15:57:58 -02:00
|
|
|
|
2012-07-28 14:11:07 -03:00
|
|
|
! set the initial density and pressure
|
2012-07-27 22:13:52 -03:00
|
|
|
!
|
2012-07-28 14:11:07 -03:00
|
|
|
q(idn,i) = dn_amb
|
2012-07-27 22:13:52 -03:00
|
|
|
#ifdef ADI
|
2012-07-28 14:11:07 -03:00
|
|
|
q(ipr,i) = pr_amb
|
2012-07-27 22:13:52 -03:00
|
|
|
#endif /* ADI */
|
2012-07-28 14:11:07 -03:00
|
|
|
|
|
|
|
! set the initial pressure in cells laying completely within the radius
|
|
|
|
!
|
|
|
|
if (ru .le. rad) then
|
|
|
|
|
|
|
|
! set the overpressure region density
|
|
|
|
!
|
|
|
|
q(idn,i) = dn_ovr
|
|
|
|
|
|
|
|
#ifdef ADI
|
|
|
|
! set the overpressure region pressure
|
|
|
|
!
|
|
|
|
q(ipr,i) = pr_ovr
|
|
|
|
#endif /* ADI */
|
|
|
|
|
|
|
|
! set the initial pressure in the cell completely outside the radius
|
|
|
|
!
|
|
|
|
else if (rl .ge. rad) then
|
|
|
|
|
|
|
|
! set the ambient region density
|
|
|
|
!
|
|
|
|
q(idn,i) = dn_amb
|
|
|
|
|
|
|
|
#ifdef ADI
|
|
|
|
! set the ambient medium pressure
|
|
|
|
!
|
|
|
|
q(ipr,i) = pr_amb
|
|
|
|
#endif /* ADI */
|
|
|
|
|
|
|
|
! integrate density or pressure in cells which are crossed by the circule with
|
|
|
|
! the given radius
|
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
else
|
2012-07-28 14:11:07 -03:00
|
|
|
|
|
|
|
#ifdef R3D
|
|
|
|
! in 3D simply set the ambient values since the integration is more complex
|
|
|
|
!
|
|
|
|
|
|
|
|
! set the ambient region density
|
|
|
|
!
|
|
|
|
q(idn,i) = dn_amb
|
|
|
|
|
2012-07-27 22:13:52 -03:00
|
|
|
#ifdef ADI
|
2012-07-28 14:11:07 -03:00
|
|
|
! set the ambient medium pressure
|
|
|
|
!
|
|
|
|
q(ipr,i) = pr_amb
|
2012-07-27 21:13:43 -03:00
|
|
|
#endif /* ADI */
|
2012-07-28 14:11:07 -03:00
|
|
|
#else /* R3D */
|
|
|
|
! calculate the bounds of area integration
|
|
|
|
!
|
|
|
|
xb = max(xl, sqrt(max(0.0d0, rad - yu * yu)))
|
|
|
|
xt = min(xu, sqrt(max(0.0d0, rad - yl * yl)))
|
|
|
|
yb = max(yl, sqrt(max(0.0d0, rad - xu * xu)))
|
|
|
|
yt = min(yu, sqrt(max(0.0d0, rad - xl * xl)))
|
|
|
|
|
|
|
|
! integrate the area below the circle within the current cell for both
|
|
|
|
! functions, y = f(x) and x = g(y), and then average them to be sure that we
|
|
|
|
! are getting the ideal symmetry
|
|
|
|
!
|
|
|
|
fc_ovr = 0.5d0 * (rad * (asin(xt / radius) - asin(xb / radius)) &
|
|
|
|
+ (xt * yb - xb * yt)) - yl * (xt - xb)
|
|
|
|
fc_ovr = fc_ovr + (xb - xl) * dy
|
|
|
|
|
|
|
|
fc_amb = 0.5d0 * (rad * (asin(yt / radius) - asin(yb / radius)) &
|
|
|
|
+ (yt * xb - yb * xt)) - xl * (yt - yb)
|
|
|
|
fc_amb = fc_amb + (yb - yl) * dx
|
|
|
|
|
|
|
|
fc_ovr = 0.5d0 * (fc_ovr + fc_amb)
|
|
|
|
|
|
|
|
! normalize coefficients
|
|
|
|
!
|
|
|
|
fc_ovr = fc_ovr / daxy
|
|
|
|
fc_amb = 1.0d0 - fc_ovr
|
|
|
|
|
|
|
|
! integrate the density over the edge cells
|
|
|
|
!
|
|
|
|
q(idn,i) = fc_ovr * dn_ovr + fc_amb * dn_amb
|
|
|
|
|
|
|
|
#ifdef ADI
|
|
|
|
! integrate the pressure over the edge cells
|
|
|
|
!
|
|
|
|
q(ipr,i) = fc_ovr * pr_ovr + fc_amb * pr_amb
|
|
|
|
#endif /* ADI */
|
|
|
|
#endif /* R3D */
|
|
|
|
|
2012-07-27 22:13:52 -03:00
|
|
|
end if
|
|
|
|
|
2008-11-29 22:19:02 -06:00
|
|
|
end do
|
|
|
|
|
2012-07-28 14:11:07 -03:00
|
|
|
! convert the primitive variables to conservative ones
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2012-07-28 14:11:07 -03:00
|
|
|
call prim2cons(im, q(1:nt,1:im), u(1:nt,1:im))
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
! copy the conserved variables to the current block
|
|
|
|
!
|
2012-07-27 22:13:52 -03:00
|
|
|
pdata%u(1:nt,1:im,j,k) = u(1:nt,1:im)
|
|
|
|
|
|
|
|
! copy the primitive variables to the current block
|
|
|
|
!
|
|
|
|
pdata%q(1:nt,1:im,j,k) = q(1:nt,1:im)
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2012-07-27 21:13:43 -03:00
|
|
|
end do
|
|
|
|
end do
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2008-12-18 10:20:15 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2012-07-27 21:13:43 -03:00
|
|
|
end subroutine setup_problem_blast
|
2008-11-11 16:12:26 -06:00
|
|
|
|
2008-12-18 10:20:15 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2012-07-27 19:23:11 -03:00
|
|
|
end module problems
|