2012-07-27 21:28:59 -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.
|
|
|
|
!!
|
2024-03-07 09:34:43 -03:00
|
|
|
!! Copyright (C) 2008-2024 Grzegorz Kowal <grzegorz@amuncode.org>
|
2012-07-27 21:28:59 -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.
|
|
|
|
!!
|
|
|
|
!! 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: REFINEMENT
|
|
|
|
!!
|
2013-12-19 14:36:37 -02:00
|
|
|
!! This module handles the error estimation and refinement criterion
|
|
|
|
!! determination.
|
2012-07-27 21:28:59 -03:00
|
|
|
!!
|
|
|
|
!!
|
|
|
|
!!******************************************************************************
|
|
|
|
!
|
|
|
|
module refinement
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
abstract interface
|
|
|
|
real(kind=4) function user_criterion_iface(pdata) result(crit)
|
|
|
|
use blocks, only : block_data
|
|
|
|
type(block_data), pointer, intent(in) :: pdata
|
|
|
|
end function
|
|
|
|
end interface
|
|
|
|
|
|
|
|
procedure(user_criterion_iface), pointer, save :: user_criterion => null()
|
|
|
|
|
2013-12-27 17:41:27 -02:00
|
|
|
real(kind=8), save :: crefmin = 2.0d-01
|
|
|
|
real(kind=8), save :: crefmax = 8.0d-01
|
2015-12-09 07:41:46 -02:00
|
|
|
real(kind=8), save :: vortmin = 1.0d-03
|
|
|
|
real(kind=8), save :: vortmax = 1.0d-01
|
2019-02-01 11:40:06 -02:00
|
|
|
real(kind=8), save :: currmin = 1.0d-03
|
|
|
|
real(kind=8), save :: currmax = 1.0d-01
|
2021-11-24 22:44:28 -03:00
|
|
|
real(kind=8), save :: usermin = 1.0d+99
|
|
|
|
real(kind=8), save :: usermax = 1.0d+99
|
2013-12-19 14:36:37 -02:00
|
|
|
real(kind=8), save :: epsref = 1.0d-02
|
|
|
|
|
2013-12-26 17:32:11 -02:00
|
|
|
logical, dimension(:), allocatable, save :: qvar_ref
|
2015-12-09 07:41:46 -02:00
|
|
|
logical , save :: vort_ref = .false.
|
2019-02-01 11:40:06 -02:00
|
|
|
logical , save :: curr_ref = .false.
|
2021-11-24 22:44:28 -03:00
|
|
|
logical , save :: user_ref = .false.
|
2012-07-27 21:28:59 -03:00
|
|
|
|
|
|
|
private
|
|
|
|
|
2019-01-30 12:47:15 -02:00
|
|
|
public :: initialize_refinement, finalize_refinement, print_refinement
|
2013-12-26 17:32:11 -02:00
|
|
|
public :: check_refinement_criterion
|
2021-11-24 22:44:28 -03:00
|
|
|
public :: user_criterion
|
2012-07-27 21:28:59 -03:00
|
|
|
|
|
|
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
|
!
|
|
|
|
contains
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!!
|
|
|
|
!!*** PUBLIC SUBROUTINES *****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-26 17:32:11 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-27 21:28:59 -03:00
|
|
|
! subroutine INITIALIZE_REFINEMENT:
|
|
|
|
! --------------------------------
|
|
|
|
!
|
2012-08-03 01:37:31 -03:00
|
|
|
! Subroutine initializes module REFINEMENT.
|
2012-07-27 21:28:59 -03:00
|
|
|
!
|
2013-12-26 17:32:11 -02:00
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! verbose - flag determining if the subroutine should be verbose;
|
2019-02-10 21:47:45 -02:00
|
|
|
! status - return flag of the procedure execution status;
|
2012-07-27 21:28:59 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2019-02-10 21:47:45 -02:00
|
|
|
subroutine initialize_refinement(verbose, status)
|
2012-07-27 21:28:59 -03:00
|
|
|
|
2019-02-10 21:47:45 -02:00
|
|
|
use equations , only : magnetized, nv, pvars
|
2021-11-24 22:44:28 -03:00
|
|
|
use helpers , only : print_message
|
2019-02-10 21:47:45 -02:00
|
|
|
use parameters, only : get_parameter
|
2012-07-27 21:28:59 -03:00
|
|
|
|
|
|
|
implicit none
|
2013-12-19 14:36:37 -02:00
|
|
|
|
2019-02-10 21:47:45 -02:00
|
|
|
logical, intent(in) :: verbose
|
|
|
|
integer, intent(out) :: status
|
2013-12-26 17:32:11 -02:00
|
|
|
|
2019-02-10 21:47:45 -02:00
|
|
|
logical :: test
|
|
|
|
integer :: p
|
|
|
|
character(len=255) :: variables = "dens pres"
|
2013-12-26 18:04:02 -02:00
|
|
|
|
2021-12-10 10:22:09 -03:00
|
|
|
character(len=*), parameter :: loc = 'REFINEMENT::initialize_refinement()'
|
2021-11-24 22:44:28 -03:00
|
|
|
|
2021-11-16 15:22:15 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2019-02-10 21:47:45 -02:00
|
|
|
!
|
|
|
|
status = 0
|
|
|
|
|
2019-01-28 21:30:45 -02:00
|
|
|
call get_parameter("crefmin", crefmin)
|
|
|
|
call get_parameter("crefmax", crefmax)
|
|
|
|
call get_parameter("vortmin", vortmin)
|
|
|
|
call get_parameter("vortmax", vortmax)
|
2019-02-01 11:40:06 -02:00
|
|
|
call get_parameter("currmin", currmin)
|
|
|
|
call get_parameter("currmax", currmax)
|
2021-11-24 22:44:28 -03:00
|
|
|
call get_parameter("usermin", usermin)
|
|
|
|
call get_parameter("usermax", usermax)
|
2019-01-28 21:30:45 -02:00
|
|
|
call get_parameter("epsref" , epsref )
|
2012-07-27 21:28:59 -03:00
|
|
|
|
2019-01-28 21:30:45 -02:00
|
|
|
call get_parameter("refinement_variables", variables)
|
2013-12-19 14:36:37 -02:00
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
allocate(qvar_ref(nv), stat=status)
|
|
|
|
if (status /= 0) &
|
|
|
|
call print_message(loc, &
|
|
|
|
"Could not allocate the space for refinement criterion names!")
|
2013-12-26 17:32:11 -02:00
|
|
|
|
|
|
|
do p = 1, nv
|
|
|
|
qvar_ref(p) = index(variables, trim(pvars(p))) > 0
|
2021-11-24 22:44:28 -03:00
|
|
|
end do
|
2013-12-26 17:32:11 -02:00
|
|
|
|
2015-12-09 07:41:46 -02:00
|
|
|
vort_ref = index(variables, 'vort') > 0
|
|
|
|
|
2019-02-10 21:47:45 -02:00
|
|
|
test = any(qvar_ref(:)) .or. vort_ref
|
|
|
|
|
2019-01-30 12:47:15 -02:00
|
|
|
if (magnetized) then
|
2019-02-01 11:40:06 -02:00
|
|
|
curr_ref = index(variables, 'curr') > 0 .or. index(variables, 'jabs') > 0
|
2019-02-10 21:47:45 -02:00
|
|
|
test = test .or. curr_ref
|
|
|
|
end if
|
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
user_ref = index(variables, 'user') > 0
|
|
|
|
test = test .or. user_ref
|
|
|
|
|
2019-02-10 21:47:45 -02:00
|
|
|
if (.not. test) then
|
2021-11-24 22:44:28 -03:00
|
|
|
if (verbose) &
|
|
|
|
call print_message(loc, "No refinement criterion has been selected!")
|
2019-02-10 21:47:45 -02:00
|
|
|
status = 1
|
2021-11-24 22:44:28 -03:00
|
|
|
return
|
2019-02-10 21:47:45 -02:00
|
|
|
end if
|
|
|
|
|
|
|
|
if (crefmin > crefmax .or. crefmin < 0.0d+00) then
|
2021-11-24 22:44:28 -03:00
|
|
|
if (verbose) &
|
|
|
|
call print_message(loc, &
|
|
|
|
"Wrong 'crefmin' or 'crefmax' parameters. " // &
|
|
|
|
"They should be positive and 'crefmin' <= 'crefmax'.")
|
2019-02-10 21:47:45 -02:00
|
|
|
status = 1
|
2021-11-24 22:44:28 -03:00
|
|
|
return
|
2019-02-10 21:47:45 -02:00
|
|
|
end if
|
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
if (vortmin > vortmax .or. vortmin < 0.0d+00) then
|
|
|
|
if (verbose) &
|
|
|
|
call print_message(loc, &
|
|
|
|
"Wrong 'vortmin' or 'vortmax' parameters. " // &
|
|
|
|
"They should be positive and 'vortmin' <= 'vortmax'.")
|
2019-02-10 21:47:45 -02:00
|
|
|
status = 1
|
2021-11-24 22:44:28 -03:00
|
|
|
return
|
2019-02-10 21:47:45 -02:00
|
|
|
end if
|
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
if ((currmin > currmax .or. currmin < 0.0d+00) .and. magnetized) then
|
|
|
|
if (verbose) &
|
|
|
|
call print_message(loc, &
|
|
|
|
"Wrong 'currmin' or 'currmax' parameters. " // &
|
|
|
|
"They should be positive and 'currmin' <= 'currmax'.")
|
2019-02-10 21:47:45 -02:00
|
|
|
status = 1
|
2021-11-24 22:44:28 -03:00
|
|
|
return
|
2019-02-10 21:47:45 -02:00
|
|
|
end if
|
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
if ((usermin > usermax .or. usermin < 0.0d+00) .and. user_ref) then
|
|
|
|
if (verbose) &
|
|
|
|
call print_message(loc, &
|
|
|
|
"Wrong 'usermin' or 'usermax' parameters. " // &
|
|
|
|
"They should be positive and 'usermin' <= 'usermax'.")
|
|
|
|
status = 1
|
|
|
|
return
|
|
|
|
end if
|
|
|
|
|
|
|
|
if (epsref <= 0.0d+00) then
|
|
|
|
if (verbose) &
|
|
|
|
call print_message(loc, &
|
|
|
|
"Wrong 'epsref' parameters. It should be positive.")
|
2019-02-10 21:47:45 -02:00
|
|
|
status = 1
|
2021-11-24 22:44:28 -03:00
|
|
|
return
|
2013-12-26 17:32:11 -02:00
|
|
|
end if
|
2013-12-19 14:36:37 -02:00
|
|
|
|
2012-07-27 21:28:59 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine initialize_refinement
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-26 17:32:11 -02:00
|
|
|
! subroutine FINALIZE_REFINEMENT:
|
|
|
|
! ------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine releases memory used by the module variables.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
2019-02-10 21:47:45 -02:00
|
|
|
! status - return flag of the procedure execution status;
|
2013-12-26 17:32:11 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2019-02-10 21:47:45 -02:00
|
|
|
subroutine finalize_refinement(status)
|
2013-12-26 17:32:11 -02:00
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
use helpers, only : print_message
|
|
|
|
|
2013-12-26 17:32:11 -02:00
|
|
|
implicit none
|
|
|
|
|
2019-02-10 21:47:45 -02:00
|
|
|
integer, intent(out) :: status
|
2013-12-26 18:04:02 -02:00
|
|
|
|
2021-12-10 10:22:09 -03:00
|
|
|
character(len=*), parameter :: loc = 'REFINEMENT:finalize_refinement()'
|
2021-11-24 22:44:28 -03:00
|
|
|
|
2021-11-16 15:22:15 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2019-02-10 21:47:45 -02:00
|
|
|
!
|
|
|
|
status = 0
|
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
if (allocated(qvar_ref)) then
|
|
|
|
deallocate(qvar_ref, stat=status)
|
|
|
|
if (status /= 0) &
|
|
|
|
call print_message(loc, &
|
|
|
|
"Could not deallocate space for the refinement criterion names!")
|
|
|
|
end if
|
2013-12-26 17:32:11 -02:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine finalize_refinement
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2019-01-30 12:47:15 -02:00
|
|
|
! subroutine PRINT_REFINEMENT:
|
|
|
|
! ---------------------------
|
|
|
|
!
|
|
|
|
! Subroutine prints module parameters.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! verbose - flag determining if the subroutine should be verbose;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine print_refinement(verbose)
|
|
|
|
|
2019-01-30 22:27:00 -02:00
|
|
|
use helpers , only : print_section, print_parameter
|
2019-01-30 12:47:15 -02:00
|
|
|
use equations, only : magnetized, pvars, nv
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
logical, intent(in) :: verbose
|
|
|
|
|
2021-11-19 11:10:38 -03:00
|
|
|
character(len=80) :: rvars = ""
|
2019-01-30 22:27:00 -02:00
|
|
|
integer :: p
|
2021-11-16 15:22:15 -03:00
|
|
|
|
2019-01-30 12:47:15 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2021-11-24 22:44:28 -03:00
|
|
|
if (.not. verbose) return
|
2019-01-30 12:47:15 -02:00
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
rvars = ""
|
|
|
|
do p = 1, nv
|
|
|
|
if (qvar_ref(p)) rvars = adjustl(trim(rvars) // ' ' // trim(pvars(p)))
|
|
|
|
end do
|
|
|
|
if (vort_ref) then
|
|
|
|
rvars = adjustl(trim(rvars) // ' vort')
|
|
|
|
end if
|
|
|
|
if (magnetized .and. curr_ref) then
|
|
|
|
rvars = adjustl(trim(rvars) // ' curr')
|
2019-01-30 12:47:15 -02:00
|
|
|
end if
|
2021-11-24 22:44:28 -03:00
|
|
|
if (user_ref) then
|
|
|
|
rvars = adjustl(trim(rvars) // ' user')
|
|
|
|
end if
|
|
|
|
|
|
|
|
call print_section(verbose, "Refinement")
|
|
|
|
call print_parameter(verbose, "refined variables", rvars)
|
|
|
|
call print_parameter(verbose, "2nd order error limits", crefmin, crefmax)
|
|
|
|
if (vort_ref) &
|
|
|
|
call print_parameter(verbose, "vorticity limits", vortmin, vortmax)
|
|
|
|
if (magnetized .and. curr_ref) &
|
|
|
|
call print_parameter(verbose, "current density limits", currmin, currmax)
|
|
|
|
if (user_ref) &
|
|
|
|
call print_parameter(verbose, "user criterion limits", usermin, usermax)
|
2019-01-30 12:47:15 -02:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine print_refinement
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-27 21:28:59 -03:00
|
|
|
! function CHECK_REFINEMENT_CRITERION:
|
|
|
|
! -----------------------------------
|
|
|
|
!
|
|
|
|
! Function scans the given data block and checks for the refinement
|
|
|
|
! criterion. It returns +1 if the criterion is met, which indicates that
|
2012-08-03 01:37:31 -03:00
|
|
|
! the block needs to be refined, 0 if there is no need for the refinement,
|
|
|
|
! and -1 if the block can be derefined.
|
2012-07-27 21:28:59 -03:00
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
2013-12-19 14:36:37 -02:00
|
|
|
! pdata - pointer to the data block for which the refinement criterion
|
|
|
|
! has to be determined;
|
2012-08-03 01:37:31 -03:00
|
|
|
!
|
2012-07-27 21:28:59 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
function check_refinement_criterion(pdata) result(criterion)
|
|
|
|
|
2021-11-16 15:22:15 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use equations, only : nv
|
2012-07-27 21:28:59 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2012-08-03 01:37:31 -03:00
|
|
|
type(block_data), pointer, intent(in) :: pdata
|
2012-07-27 21:28:59 -03:00
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
integer(kind=4) :: criterion
|
2012-07-27 21:28:59 -03:00
|
|
|
|
2013-12-26 17:32:11 -02:00
|
|
|
integer :: p
|
2020-08-06 21:39:10 -03:00
|
|
|
real(kind=8) :: cref
|
2013-12-26 18:04:02 -02:00
|
|
|
|
2021-11-16 15:22:15 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2012-07-27 21:28:59 -03:00
|
|
|
!
|
2015-12-09 07:41:46 -02:00
|
|
|
criterion = -1
|
2012-07-27 21:28:59 -03:00
|
|
|
|
2013-12-26 17:32:11 -02:00
|
|
|
do p = 1, nv
|
2015-12-09 07:41:46 -02:00
|
|
|
if (qvar_ref(p)) then
|
|
|
|
|
|
|
|
cref = second_derivative_error(p, pdata)
|
|
|
|
|
|
|
|
if (cref > crefmin) criterion = max(criterion, 0)
|
|
|
|
if (cref > crefmax) criterion = max(criterion, 1)
|
|
|
|
|
|
|
|
end if
|
2021-11-24 22:44:28 -03:00
|
|
|
end do
|
2012-07-27 21:28:59 -03:00
|
|
|
|
2015-12-09 07:41:46 -02:00
|
|
|
if (vort_ref) then
|
|
|
|
cref = vorticity_magnitude(pdata)
|
|
|
|
|
|
|
|
if (cref > vortmin) criterion = max(criterion, 0)
|
|
|
|
if (cref > vortmax) criterion = max(criterion, 1)
|
|
|
|
end if
|
2014-09-23 11:24:52 -03:00
|
|
|
|
2019-02-01 11:40:06 -02:00
|
|
|
if (curr_ref) then
|
2015-12-09 07:41:46 -02:00
|
|
|
cref = current_density_magnitude(pdata)
|
2012-07-27 21:28:59 -03:00
|
|
|
|
2019-02-01 11:40:06 -02:00
|
|
|
if (cref > currmin) criterion = max(criterion, 0)
|
|
|
|
if (cref > currmax) criterion = max(criterion, 1)
|
2012-07-27 21:28:59 -03:00
|
|
|
end if
|
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
if (user_ref .and. associated(user_criterion)) then
|
|
|
|
cref = user_criterion(pdata)
|
|
|
|
|
|
|
|
if (cref > usermin) criterion = max(criterion, 0)
|
|
|
|
if (cref > usermax) criterion = max(criterion, 1)
|
|
|
|
end if
|
|
|
|
|
2012-07-27 21:28:59 -03:00
|
|
|
return
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function check_refinement_criterion
|
2013-12-19 14:36:37 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!!
|
|
|
|
!!*** PRIVATE SUBROUTINES ****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-26 17:32:11 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-19 14:36:37 -02:00
|
|
|
! function SECOND_DERIVATIVE_ERROR:
|
|
|
|
! --------------------------------
|
|
|
|
!
|
|
|
|
! Function calculate the second derivative error for a given data block
|
|
|
|
! and selected primitive variables. The total error is returned then.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! iqt - the index of primitive variable;
|
|
|
|
! pdata - pointer to the data block for which error is calculated;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
function second_derivative_error(iqt, pdata) result(error)
|
|
|
|
|
2019-02-05 12:59:44 -02:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : nbl, neu
|
2013-12-19 14:36:37 -02:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer , intent(in) :: iqt
|
|
|
|
type(block_data), pointer, intent(in) :: pdata
|
|
|
|
|
2020-08-06 21:39:10 -03:00
|
|
|
real(kind=8) :: error
|
2013-12-19 14:36:37 -02:00
|
|
|
|
|
|
|
integer :: i, im1, ip1
|
|
|
|
integer :: j, jm1, jp1
|
2022-01-08 11:02:05 -03:00
|
|
|
integer :: k
|
2020-08-15 01:27:08 -03:00
|
|
|
#if NDIMS == 3
|
|
|
|
integer :: km1, kp1
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
real(kind=8) :: fl, fr, fc, fx, fy
|
|
|
|
#if NDIMS == 3
|
|
|
|
real(kind=8) :: fz
|
|
|
|
#endif /* NDIMS == 3 */
|
2013-12-19 14:51:58 -02:00
|
|
|
|
|
|
|
real(kind=8), parameter :: eps = epsilon(1.0d+00)
|
2021-11-16 15:22:15 -03:00
|
|
|
|
2013-12-19 14:36:37 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2014-01-31 10:50:18 -02:00
|
|
|
error = 0.0e+00
|
2013-12-19 14:36:37 -02:00
|
|
|
|
2022-01-08 11:02:05 -03:00
|
|
|
#if NDIMS == 2
|
|
|
|
k = 1
|
|
|
|
#endif /* NDIMS == 2 */
|
|
|
|
|
2013-12-19 14:36:37 -02:00
|
|
|
if (iqt > 0) then
|
|
|
|
|
|
|
|
#if NDIMS == 3
|
2019-02-05 12:59:44 -02:00
|
|
|
do k = nbl, neu
|
2013-12-19 14:36:37 -02:00
|
|
|
km1 = k - 1
|
|
|
|
kp1 = k + 1
|
|
|
|
#endif /* NDIMS == 3 */
|
2019-02-05 12:59:44 -02:00
|
|
|
do j = nbl, neu
|
2013-12-19 14:36:37 -02:00
|
|
|
jm1 = j - 1
|
|
|
|
jp1 = j + 1
|
2019-02-05 12:59:44 -02:00
|
|
|
do i = nbl, neu
|
2013-12-19 14:36:37 -02:00
|
|
|
im1 = i - 1
|
|
|
|
ip1 = i + 1
|
|
|
|
|
|
|
|
fr = pdata%q(iqt,ip1,j,k) - pdata%q(iqt,i ,j,k)
|
|
|
|
fl = pdata%q(iqt,im1,j,k) - pdata%q(iqt,i ,j,k)
|
2021-11-24 22:44:28 -03:00
|
|
|
fc = abs(pdata%q(iqt,ip1,j,k)) + abs(pdata%q(iqt,im1,j,k)) &
|
|
|
|
+ 2.0d+00 * abs(pdata%q(iqt,i,j,k))
|
2014-01-31 10:50:18 -02:00
|
|
|
fx = abs(fr + fl) / (abs(fr) + abs(fl) + epsref * fc + eps)
|
2013-12-19 14:36:37 -02:00
|
|
|
|
|
|
|
fr = pdata%q(iqt,i,jp1,k) - pdata%q(iqt,i,j ,k)
|
|
|
|
fl = pdata%q(iqt,i,jm1,k) - pdata%q(iqt,i,j ,k)
|
2021-11-24 22:44:28 -03:00
|
|
|
fc = abs(pdata%q(iqt,i,jp1,k)) + abs(pdata%q(iqt,i,jm1,k)) &
|
|
|
|
+ 2.0d+00 * abs(pdata%q(iqt,i,j,k))
|
2014-01-31 10:50:18 -02:00
|
|
|
fy = abs(fr + fl) / (abs(fr) + abs(fl) + epsref * fc + eps)
|
2013-12-19 14:36:37 -02:00
|
|
|
|
|
|
|
#if NDIMS == 3
|
|
|
|
fr = pdata%q(iqt,i,j,kp1) - pdata%q(iqt,i,j,k )
|
|
|
|
fl = pdata%q(iqt,i,j,km1) - pdata%q(iqt,i,j,k )
|
2021-11-24 22:44:28 -03:00
|
|
|
fc = abs(pdata%q(iqt,i,j,kp1)) + abs(pdata%q(iqt,i,j,km1)) &
|
|
|
|
+ 2.0d+00 * abs(pdata%q(iqt,i,j,k))
|
2014-01-31 10:50:18 -02:00
|
|
|
fz = abs(fr + fl) / (abs(fr) + abs(fl) + epsref * fc + eps)
|
2013-12-19 14:36:37 -02:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
#if NDIMS == 2
|
2014-01-31 10:50:18 -02:00
|
|
|
error = max(error, fx, fy)
|
2013-12-19 14:36:37 -02:00
|
|
|
#endif /* NDIMS == 2 */
|
|
|
|
#if NDIMS == 3
|
2014-01-31 10:50:18 -02:00
|
|
|
error = max(error, fx, fy, fz)
|
2013-12-19 14:36:37 -02:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
end do
|
|
|
|
end do
|
2019-02-05 12:59:44 -02:00
|
|
|
#if NDIMS == 3
|
2021-11-24 22:44:28 -03:00
|
|
|
end do
|
2019-02-05 12:59:44 -02:00
|
|
|
#endif /* NDIMS == 3 */
|
2013-12-19 14:36:37 -02:00
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
end if
|
2013-12-19 14:36:37 -02:00
|
|
|
|
|
|
|
return
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function second_derivative_error
|
2014-09-23 11:24:52 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2015-12-09 07:41:46 -02:00
|
|
|
! function VORTICITY_MAGNITUDE:
|
|
|
|
! ----------------------------
|
|
|
|
!
|
2021-11-11 17:54:54 -03:00
|
|
|
! Function finds the maximum values of the vorticity magnitude
|
|
|
|
! for the current data block.
|
2015-12-09 07:41:46 -02:00
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pdata - pointer to the data block for which error is calculated;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
function vorticity_magnitude(pdata) result(wmax)
|
|
|
|
|
2021-11-19 12:44:56 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : nn => bcells
|
|
|
|
use coordinates, only : nbl, neu
|
|
|
|
use equations , only : ivx, ivz
|
|
|
|
use helpers , only : print_message
|
|
|
|
use operators , only : curl
|
|
|
|
use workspace , only : resize_workspace, work, work_in_use
|
2015-12-09 07:41:46 -02:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
type(block_data), pointer, intent(in) :: pdata
|
|
|
|
|
|
|
|
real(kind=4) :: wmax
|
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
logical, save :: first = .true.
|
|
|
|
|
2022-01-08 11:02:05 -03:00
|
|
|
integer :: i, j, k, status
|
2015-12-09 07:41:46 -02:00
|
|
|
real(kind=8) :: vort
|
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
real(kind=8), dimension(3), save :: dh
|
|
|
|
|
|
|
|
real(kind=8), dimension(:,:,:,:), pointer, save :: w
|
|
|
|
|
2022-01-08 11:52:16 -03:00
|
|
|
integer, save :: nt
|
2021-12-07 10:46:18 -03:00
|
|
|
!$ integer :: omp_get_thread_num
|
2021-12-10 10:22:09 -03:00
|
|
|
!$omp threadprivate(first,nt,dh,w)
|
2021-12-07 10:46:18 -03:00
|
|
|
|
2021-12-10 10:22:09 -03:00
|
|
|
character(len=*), parameter :: loc = 'REFINEMENT::vorticity_magnitude()'
|
2021-11-11 17:54:54 -03:00
|
|
|
|
2015-12-09 07:41:46 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2022-01-08 11:52:16 -03:00
|
|
|
nt = 0
|
2021-12-07 19:55:30 -03:00
|
|
|
!$ nt = omp_get_thread_num()
|
2021-11-12 14:08:04 -03:00
|
|
|
wmax = 0.0e+00
|
2021-11-12 13:55:29 -03:00
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
if (first) then
|
2021-11-13 20:29:08 -03:00
|
|
|
i = 3 * nn**NDIMS
|
|
|
|
|
|
|
|
call resize_workspace(i, status)
|
|
|
|
if (status /= 0) then
|
2021-11-19 12:44:56 -03:00
|
|
|
call print_message(loc, "Could not resize the workspace!")
|
2021-11-12 22:42:53 -03:00
|
|
|
return
|
2021-11-13 20:29:08 -03:00
|
|
|
end if
|
2021-11-12 22:42:53 -03:00
|
|
|
|
|
|
|
#if NDIMS == 3
|
2021-12-07 10:46:18 -03:00
|
|
|
w(1:3,1:nn,1:nn,1:nn) => work(1:i,nt)
|
2021-11-12 22:42:53 -03:00
|
|
|
#else /* NDIMS == 3 */
|
2021-12-07 10:46:18 -03:00
|
|
|
w(1:3,1:nn,1:nn,1: 1) => work(1:i,nt)
|
2021-11-12 22:42:53 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
2021-11-13 20:29:08 -03:00
|
|
|
dh(:) = 1.0d+00
|
2021-11-12 22:42:53 -03:00
|
|
|
|
|
|
|
first = .false.
|
|
|
|
end if
|
|
|
|
|
2021-12-07 10:46:18 -03:00
|
|
|
if (work_in_use(nt)) &
|
2021-11-24 22:44:28 -03:00
|
|
|
call print_message(loc, &
|
|
|
|
"Workspace is being used right now! Corruptions can occur!")
|
2021-11-13 20:29:08 -03:00
|
|
|
|
2021-12-07 10:46:18 -03:00
|
|
|
work_in_use(nt) = .true.
|
2021-11-12 22:42:53 -03:00
|
|
|
|
|
|
|
call curl(dh(:), pdata%q(ivx:ivz,:,:,:), w(:,:,:,:))
|
|
|
|
|
2022-01-08 11:02:05 -03:00
|
|
|
#if NDIMS == 2
|
|
|
|
k = 1
|
|
|
|
#endif /* NDIMS == 2 */
|
2019-02-05 11:06:01 -02:00
|
|
|
#if NDIMS == 3
|
2021-11-11 17:54:54 -03:00
|
|
|
do k = nbl, neu
|
2019-02-05 11:06:01 -02:00
|
|
|
#endif /* NDIMS == 3 */
|
2021-11-11 17:54:54 -03:00
|
|
|
do j = nbl, neu
|
|
|
|
do i = nbl, neu
|
2015-12-09 07:41:46 -02:00
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
vort = sum(w(:,i,j,k)**2)
|
2015-12-09 07:41:46 -02:00
|
|
|
|
2021-11-11 17:54:54 -03:00
|
|
|
wmax = max(wmax, real(vort, kind=4))
|
2015-12-09 07:41:46 -02:00
|
|
|
|
2021-11-11 17:54:54 -03:00
|
|
|
end do
|
|
|
|
end do
|
2019-02-05 11:06:01 -02:00
|
|
|
#if NDIMS == 3
|
2021-11-11 17:54:54 -03:00
|
|
|
end do
|
2019-02-05 11:06:01 -02:00
|
|
|
#endif /* NDIMS == 3 */
|
2015-12-09 07:41:46 -02:00
|
|
|
|
2021-12-07 10:46:18 -03:00
|
|
|
work_in_use(nt) = .false.
|
2021-11-12 22:42:53 -03:00
|
|
|
|
2021-11-11 17:54:54 -03:00
|
|
|
wmax = sqrt(wmax)
|
2015-12-09 07:41:46 -02:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function vorticity_magnitude
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2014-09-23 11:24:52 -03:00
|
|
|
! function CURRENT_DENSITY_MAGNITUDE:
|
|
|
|
! ----------------------------------
|
|
|
|
!
|
2021-11-11 17:54:54 -03:00
|
|
|
! Function finds the maximum values of the current density magnitude
|
|
|
|
! for the current data block.
|
2014-09-23 11:24:52 -03:00
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pdata - pointer to the data block for which error is calculated;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
function current_density_magnitude(pdata) result(jmax)
|
|
|
|
|
2021-11-19 12:44:56 -03:00
|
|
|
use blocks , only : block_data
|
|
|
|
use coordinates, only : nn => bcells
|
|
|
|
use coordinates, only : nbl, neu
|
|
|
|
use equations , only : magnetized
|
|
|
|
use equations , only : ibx, ibz
|
|
|
|
use helpers , only : print_message
|
|
|
|
use operators , only : curl
|
|
|
|
use workspace , only : resize_workspace, work, work_in_use
|
2014-09-23 11:24:52 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
type(block_data), pointer, intent(in) :: pdata
|
|
|
|
|
|
|
|
real(kind=4) :: jmax
|
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
logical, save :: first = .true.
|
|
|
|
|
2022-01-08 11:02:05 -03:00
|
|
|
integer :: i, j, k, status
|
2014-09-23 11:24:52 -03:00
|
|
|
real(kind=8) :: jabs
|
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
real(kind=8), dimension(3), save :: dh
|
|
|
|
|
|
|
|
real(kind=8), dimension(:,:,:,:), pointer, save :: w
|
|
|
|
|
2022-01-08 11:52:16 -03:00
|
|
|
integer, save :: nt
|
2021-12-07 10:46:18 -03:00
|
|
|
!$ integer :: omp_get_thread_num
|
2021-12-10 10:22:09 -03:00
|
|
|
!$omp threadprivate(first,nt,dh,w)
|
2021-12-07 10:46:18 -03:00
|
|
|
|
2021-12-10 10:22:09 -03:00
|
|
|
character(len=*), parameter :: loc = &
|
|
|
|
'REFINEMENT::current_density_magnitude()'
|
2021-11-11 17:54:54 -03:00
|
|
|
|
2014-09-23 11:24:52 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
jmax = 0.0e+00
|
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
if (.not. magnetized) return
|
2022-01-08 11:52:16 -03:00
|
|
|
|
|
|
|
nt = 0
|
2021-12-10 10:22:09 -03:00
|
|
|
!$ nt = omp_get_thread_num()
|
2021-11-12 22:42:53 -03:00
|
|
|
|
|
|
|
if (first) then
|
2021-11-13 19:46:08 -03:00
|
|
|
i = 3 * nn**NDIMS
|
2021-11-12 22:42:53 -03:00
|
|
|
|
2021-11-13 19:46:08 -03:00
|
|
|
call resize_workspace(i, status)
|
|
|
|
if (status /= 0) then
|
2021-11-19 12:44:56 -03:00
|
|
|
call print_message(loc, "Could not resize the workspace!")
|
2021-11-13 19:46:08 -03:00
|
|
|
return
|
|
|
|
end if
|
2021-11-12 22:42:53 -03:00
|
|
|
|
|
|
|
#if NDIMS == 3
|
2021-12-07 10:46:18 -03:00
|
|
|
w(1:3,1:nn,1:nn,1:nn) => work(1:i,nt)
|
2021-11-12 22:42:53 -03:00
|
|
|
#else /* NDIMS == 3 */
|
2021-12-07 10:46:18 -03:00
|
|
|
w(1:3,1:nn,1:nn,1: 1) => work(1:i,nt)
|
2021-11-12 22:42:53 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
2021-11-13 19:46:08 -03:00
|
|
|
dh(:) = 1.0d+00
|
2021-11-12 22:42:53 -03:00
|
|
|
|
|
|
|
first = .false.
|
|
|
|
end if
|
|
|
|
|
2021-12-07 10:46:18 -03:00
|
|
|
if (work_in_use(nt)) &
|
2021-11-24 22:44:28 -03:00
|
|
|
call print_message(loc, &
|
|
|
|
"Workspace is being used right now! Corruptions can occur!")
|
2021-11-13 20:26:34 -03:00
|
|
|
|
2021-12-07 10:46:18 -03:00
|
|
|
work_in_use(nt) = .true.
|
2014-09-23 11:24:52 -03:00
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
call curl(dh(:), pdata%q(ibx:ibz,:,:,:), w(:,:,:,:))
|
2014-09-23 11:24:52 -03:00
|
|
|
|
2022-01-08 11:02:05 -03:00
|
|
|
#if NDIMS == 2
|
|
|
|
k = 1
|
|
|
|
#endif /* NDIMS == 2 */
|
2019-02-05 11:06:01 -02:00
|
|
|
#if NDIMS == 3
|
2021-11-12 22:42:53 -03:00
|
|
|
do k = nbl, neu
|
2019-02-05 11:06:01 -02:00
|
|
|
#endif /* NDIMS == 3 */
|
2021-11-12 22:42:53 -03:00
|
|
|
do j = nbl, neu
|
|
|
|
do i = nbl, neu
|
2014-09-23 11:24:52 -03:00
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
jabs = sum(w(:,i,j,k)**2)
|
2014-09-23 11:24:52 -03:00
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
jmax = max(jmax, real(jabs, kind=4))
|
2014-09-23 11:24:52 -03:00
|
|
|
|
2021-11-24 22:44:28 -03:00
|
|
|
end do
|
|
|
|
end do
|
2019-02-05 11:06:01 -02:00
|
|
|
#if NDIMS == 3
|
2021-11-24 22:44:28 -03:00
|
|
|
end do
|
2019-02-05 11:06:01 -02:00
|
|
|
#endif /* NDIMS == 3 */
|
2014-09-23 11:24:52 -03:00
|
|
|
|
2021-12-07 10:46:18 -03:00
|
|
|
work_in_use(nt) = .false.
|
2014-09-23 14:17:38 -03:00
|
|
|
|
2021-11-12 22:42:53 -03:00
|
|
|
jmax = sqrt(jmax)
|
2021-11-11 17:54:54 -03:00
|
|
|
|
2014-09-23 11:24:52 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function current_density_magnitude
|
2012-07-27 21:28:59 -03:00
|
|
|
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
end module refinement
|