amun-code/sources/user_problem.F90
2021-11-18 13:08:01 -03:00

462 lines
14 KiB
Fortran

!!******************************************************************************
!!
!! This file is part of the AMUN source code, a program to perform
!! Newtonian or relativistic magnetohydrodynamical simulations on uniform or
!! adaptive mesh.
!!
!! Copyright (C) 2017-2021 Grzegorz Kowal <grzegorz@amuncode.org>
!!
!! 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.
!!
!! This program is distributed in the hope that it will be useful,
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
!! GNU General Public License for more details.
!!
!! You should have received a copy of the GNU General Public License
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
!!
!!*****************************************************************************
!!
!! module: USER_PROBLEM
!!
!! This module provides subroutines to setup custom problem.
!!
!!*****************************************************************************
!
module user_problem
implicit none
! default problem parameters
!
real(kind=8), save :: dens = 1.0d+00
real(kind=8), save :: beta = 1.0d+00
real(kind=8), save :: pres = 1.0d+00
real(kind=8), save :: bmag = sqrt(2.0d+00)
! flag indicating if the gravitational source term should be enabled
!
logical, save :: user_gravity_enabled = .false.
private
public :: initialize_user_problem, finalize_user_problem
public :: setup_user_problem, update_user_sources, update_user_shapes
public :: user_gravitational_acceleration, user_gravity_enabled
public :: user_boundary_x, user_boundary_y, user_boundary_z
public :: user_statistics
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!
contains
!
!===============================================================================
!
! subroutine INITIALIZE_USER_PROBLEM:
! ----------------------------------
!
! Subroutine initializes user problem. It could read problem parameters which
! are used in all subroutines defining this specific problem.
!
! Arguments:
!
! problem - the problem name
! rcount - the run count for restarted jobs
! verbose - a logical flag turning the information printing;
! status - an integer flag for error return value;
!
!===============================================================================
!
subroutine initialize_user_problem(problem, rcount, verbose, status)
use equations , only : eos, csnd, csnd2, adiabatic_index
use helpers , only : print_section, print_parameter
use parameters, only : get_parameter
implicit none
character(len=64), intent(in) :: problem
integer , intent(in) :: rcount
logical , intent(in) :: verbose
integer , intent(out) :: status
!-------------------------------------------------------------------------------
!
status = 0
! get parameters
!
call get_parameter("density", dens)
call get_parameter("beta" , beta)
if (eos == "iso") then
pres = csnd2 * dens
else
pres = csnd2 * dens / adiabatic_index
end if
bmag = sqrt(2.0d+00 * pres / beta)
! print information about the user problem setup
!
call print_section(verbose, "Parameters (* - set, otherwise calculated)")
call print_parameter(verbose, '(*) beta (plasma-beta)', beta)
call print_parameter(verbose, '(*) dens' , dens)
call print_parameter(verbose, '( ) pres' , pres)
call print_parameter(verbose, '( ) bmag' , bmag)
call print_parameter(verbose, '(*) csnd' , csnd)
call print_parameter(verbose, '( ) valf' , bmag / sqrt(dens))
!-------------------------------------------------------------------------------
!
end subroutine initialize_user_problem
!
!===============================================================================
!
! subroutine FINALIZE_USER_PROBLEM:
! --------------------------------
!
! Subroutine releases memory used by the module.
!
! Arguments:
!
! status - an integer flag for error return value;
!
!===============================================================================
!
subroutine finalize_user_problem(status)
implicit none
integer, intent(out) :: status
!-------------------------------------------------------------------------------
!
status = 0
!-------------------------------------------------------------------------------
!
end subroutine finalize_user_problem
!
!===============================================================================
!
! subroutine SETUP_USER_PROBLEM:
! -----------------------------
!
! Subroutine sets the initial conditions for the user specific problem.
!
! Arguments:
!
! pdata - pointer to the datablock structure of the currently initialized
! block;
!
!===============================================================================
!
subroutine setup_user_problem(pdata)
use blocks , only : block_data
use coordinates, only : nn => bcells
use equations , only : prim2cons
use equations , only : nv
use equations , only : idn, ivx, ivy, ivz, ipr, ibx, iby, ibz, ibp
implicit none
type(block_data), pointer, intent(inout) :: pdata
integer :: j, k = 1
real(kind=8), dimension(nv,nn) :: q, u
!-------------------------------------------------------------------------------
!
! set the variables
!
q(idn,:) = dens
if (ipr > 0) q(ipr,:) = pres
q(ivx,:) = 0.0d+00
q(ivy,:) = 0.0d+00
q(ivz,:) = 0.0d+00
if (ibx > 0) then
q(ibx,:) = bmag
q(iby,:) = 0.0d+00
q(ibz,:) = 0.0d+00
q(ibp,:) = 0.0d+00
end if
! convert the primitive variables to conservative ones
!
call prim2cons(q(:,:), u(:,:))
! iterate over all positions in the YZ plane
!
#if NDIMS == 3
do k = 1, nn
#endif /* NDIMS == 3 */
do j = 1, nn
! copy the primitive variables to the current block
!
pdata%q(:,:,j,k) = q(:,:)
! copy the conserved variables to the current block
!
pdata%u(:,:,j,k) = u(:,:)
end do ! j = 1, nn
#if NDIMS == 3
end do ! k = 1, nn
#endif /* NDIMS == 3 */
!-------------------------------------------------------------------------------
!
end subroutine setup_user_problem
!
!===============================================================================
!
! subroutine UPDATE_USER_SOURCES:
! ------------------------------
!
! Subroutine adds the user defined source terms.
!
! Arguments:
!
! pdata - the pointer to a data block;
! t, dt - the time and time increment;
! du - the array of variable increment;
!
!===============================================================================
!
subroutine update_user_sources(pdata, t, dt, du)
use blocks, only : block_data
implicit none
type(block_data), pointer , intent(inout) :: pdata
real(kind=8) , intent(in) :: t, dt
real(kind=8), dimension(:,:,:,:), intent(inout) :: du
!-------------------------------------------------------------------------------
!
!-------------------------------------------------------------------------------
!
end subroutine update_user_sources
!
!===============================================================================
!
! subroutine UPDATE_USER_SHAPES:
! -----------------------------
!
! Subroutine defines the regions updated by user.
!
! Arguments:
!
! pdata - pointer to the data block structure of the currently initialized
! block;
! time - time at the moment of update;
! dt - time step since the last update;
!
!===============================================================================
!
subroutine update_user_shapes(pdata, time, dt)
use blocks, only : block_data
implicit none
type(block_data), pointer, intent(inout) :: pdata
real(kind=8) , intent(in) :: time, dt
!-------------------------------------------------------------------------------
!
!-------------------------------------------------------------------------------
!
end subroutine update_user_shapes
!
!===============================================================================
!
! subroutine USER_GRAVITATIONAL_ACCELERATION:
! ------------------------------------------
!
! Subroutine returns the user defined gravitational acceleration.
!
! Arguments:
!
! t, dt - time and the time increment;
! x, y, z - rectangular coordinates;
! acc - vector of the gravitational acceleration;
!
!===============================================================================
!
subroutine user_gravitational_acceleration(t, dt, x, y, z, acc)
use parameters, only : get_parameter
implicit none
real(kind=8) , intent(in) :: t, dt
real(kind=8) , intent(in) :: x, y, z
real(kind=8), dimension(3), intent(out) :: acc
!-------------------------------------------------------------------------------
!
acc(:) = 0.0d+00
!-------------------------------------------------------------------------------
!
end subroutine user_gravitational_acceleration
!
!===============================================================================
!
! subroutine USER_BOUNDARY_X:
! --------------------------
!
! Subroutine updates ghost zones within the specific region along
! the X direction.
!
! Arguments:
!
! ic - the block side along the X direction for the ghost zone update;
! jl, ju - the cell index limits for the Y direction;
! kl, ku - the cell index limits for the Z direction;
! t, dt - time and time increment;
! x, y, z - the block coordinates;
! qn - the array of variables to update;
!
!===============================================================================
!
subroutine user_boundary_x(ic, jl, ju, kl, ku, t, dt, x, y, z, qn)
implicit none
integer , intent(in) :: ic
integer , intent(in) :: jl, ju
integer , intent(in) :: kl, ku
real(kind=8) , intent(in) :: t, dt
real(kind=8), dimension(:) , intent(in) :: x
real(kind=8), dimension(:) , intent(in) :: y
real(kind=8), dimension(:) , intent(in) :: z
real(kind=8), dimension(:,:,:,:), intent(inout) :: qn
!-------------------------------------------------------------------------------
!
!-------------------------------------------------------------------------------
!
end subroutine user_boundary_x
!
!===============================================================================
!
! subroutine USER_BOUNDARY_Y:
! --------------------------
!
! Subroutine updates ghost zones within the specific region along
! the Y direction.
!
! Arguments:
!
! jc - the block side along the Y direction for the ghost zone update;
! il, iu - the cell index limits for the X direction;
! kl, ku - the cell index limits for the Z direction;
! t, dt - time and time increment;
! x, y, z - the block coordinates;
! qn - the array of variables to update;
!
!===============================================================================
!
subroutine user_boundary_y(jc, il, iu, kl, ku, t, dt, x, y, z, qn)
implicit none
integer , intent(in) :: jc
integer , intent(in) :: il, iu
integer , intent(in) :: kl, ku
real(kind=8) , intent(in) :: t, dt
real(kind=8), dimension(:) , intent(in) :: x
real(kind=8), dimension(:) , intent(in) :: y
real(kind=8), dimension(:) , intent(in) :: z
real(kind=8), dimension(:,:,:,:), intent(inout) :: qn
!-------------------------------------------------------------------------------
!
!-------------------------------------------------------------------------------
!
end subroutine user_boundary_y
!
!===============================================================================
!
! subroutine USER_BOUNDARY_Z:
! --------------------------
!
! Subroutine updates ghost zones within the specific region along
! the Z direction.
!
! Arguments:
!
! kc - the block side along the Z direction for the ghost zone update;
! il, iu - the cell index limits for the X direction;
! jl, ju - the cell index limits for the Y direction;
! t, dt - time and time increment;
! x, y, z - the block coordinates;
! qn - the array of variables to update;
!
!===============================================================================
!
subroutine user_boundary_z(kc, il, iu, jl, ju, t, dt, x, y, z, qn)
implicit none
integer , intent(in) :: kc
integer , intent(in) :: il, iu
integer , intent(in) :: jl, ju
real(kind=8) , intent(in) :: t, dt
real(kind=8), dimension(:) , intent(in) :: x
real(kind=8), dimension(:) , intent(in) :: y
real(kind=8), dimension(:) , intent(in) :: z
real(kind=8), dimension(:,:,:,:), intent(inout) :: qn
!-------------------------------------------------------------------------------
!
!-------------------------------------------------------------------------------
!
end subroutine user_boundary_z
!
!===============================================================================
!
! subroutine USER_STATISTICS:
! -------------------------------
!
! Subroutine can be use to store user defined time statistics. The file to
! store these statistics should be properly created in subroutine
! initialize_user_problem() and closed in finalize_user_problem().
!
! Arguments:
!
! time - the current simulation time;
!
!===============================================================================
!
subroutine user_statistics(time)
implicit none
real(kind=8), intent(in) :: time
!-------------------------------------------------------------------------------
!
!-------------------------------------------------------------------------------
!
end subroutine user_statistics
!===============================================================================
!
end module user_problem