2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2008-11-11 16:12:26 -06:00
|
|
|
!! module: blocks - handling block storage
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2011-02-27 22:45:54 -03:00
|
|
|
!! Copyright (C) 2008-2011 Grzegorz Kowal <grzegorz@gkowal.info>
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2011-04-25 13:44:34 -03:00
|
|
|
!! This file is part of AMUN.
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2011-04-29 11:21:30 -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 2
|
|
|
|
!! of the License, or (at your option) any later version.
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2011-04-29 11:21:30 -03:00
|
|
|
!! This program is distributed in the hope that it will be useful,
|
2008-11-04 21:00:50 -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
|
2011-04-29 11:21:30 -03:00
|
|
|
!! along with this program; if not, write to the Free Software Foundation,
|
|
|
|
!! Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
|
|
|
!
|
|
|
|
module blocks
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! parameters
|
|
|
|
!
|
2008-12-08 15:31:35 -06:00
|
|
|
integer(kind=4), parameter :: ndims = NDIMS
|
2009-09-11 21:52:18 -03:00
|
|
|
integer(kind=4), parameter :: nsides = 2
|
|
|
|
integer(kind=4), parameter :: nfaces = 2**(ndims-1)
|
2008-11-04 21:00:50 -06:00
|
|
|
integer(kind=4), parameter :: nchild = 2**ndims
|
|
|
|
|
2009-09-18 17:43:48 -03:00
|
|
|
!! BLOCK STRUCTURE POINTERS (have to be defined before structures)
|
|
|
|
!!
|
2009-09-09 16:37:28 -03:00
|
|
|
! define pointers to block_meta and block_data structures
|
|
|
|
!
|
|
|
|
type pointer_meta
|
|
|
|
type(block_meta), pointer :: ptr
|
|
|
|
end type pointer_meta
|
|
|
|
|
|
|
|
type pointer_data
|
|
|
|
type(block_data), pointer :: ptr
|
|
|
|
end type pointer_data
|
|
|
|
|
2009-09-18 17:43:48 -03:00
|
|
|
type pointer_info
|
|
|
|
type(block_info), pointer :: ptr
|
|
|
|
end type pointer_info
|
|
|
|
|
|
|
|
!! BLOCK STRUCTURES
|
|
|
|
!
|
2009-09-09 16:37:28 -03:00
|
|
|
! define block_meta structure
|
|
|
|
!
|
|
|
|
type block_meta
|
|
|
|
type(block_meta) , pointer :: prev ! pointer to the previous block
|
|
|
|
type(block_meta) , pointer :: next ! pointer to the next block
|
|
|
|
type(block_meta) , pointer :: parent ! pointer to the parent block
|
|
|
|
type(pointer_meta) :: child(nchild) ! pointers to children
|
2009-09-26 14:27:47 -03:00
|
|
|
type(pointer_meta) :: neigh(ndims,nsides,nfaces)
|
|
|
|
! pointers to neighbors
|
2009-09-09 16:37:28 -03:00
|
|
|
|
|
|
|
type(block_data) , pointer :: data ! pointer to the data block
|
|
|
|
|
|
|
|
integer(kind=4) :: id ! block identificator
|
|
|
|
integer(kind=4) :: cpu ! the cpu id of the block
|
|
|
|
integer(kind=4) :: level ! refinement level
|
2009-09-10 17:25:28 -03:00
|
|
|
integer(kind=4) :: config ! configuration flag
|
2009-09-09 16:37:28 -03:00
|
|
|
integer(kind=4) :: refine ! refinement flag:
|
|
|
|
! -1 - derefine
|
|
|
|
! 0 - do nothing
|
|
|
|
! 1 - refine
|
|
|
|
|
2010-12-01 09:31:38 -02:00
|
|
|
integer(kind=4) :: pos(ndims) ! the position in the parent block
|
2010-10-11 02:10:47 -03:00
|
|
|
integer(kind=4) :: coord(ndims) ! coordinate of the lower
|
|
|
|
! corner in the effective
|
|
|
|
! resolution
|
|
|
|
|
2009-09-09 16:37:28 -03:00
|
|
|
logical :: leaf ! leaf flag
|
2009-09-25 11:41:19 -03:00
|
|
|
|
2010-10-11 01:39:15 -03:00
|
|
|
real(kind=8) :: xmin, xmax ! bounds for the x direction
|
|
|
|
real(kind=8) :: ymin, ymax ! bounds for the y direction
|
|
|
|
real(kind=8) :: zmin, zmax ! bounds for the z direction
|
2009-09-09 16:37:28 -03:00
|
|
|
end type block_meta
|
|
|
|
|
|
|
|
! define block_data structure
|
|
|
|
!
|
|
|
|
type block_data
|
|
|
|
type(block_data), pointer :: prev ! pointer to the previous block
|
|
|
|
type(block_data), pointer :: next ! pointer to the next block
|
|
|
|
|
|
|
|
type(block_meta), pointer :: meta ! pointer to the metadata block
|
|
|
|
|
2011-04-30 12:11:57 -03:00
|
|
|
real, dimension(:,:,:,:) , allocatable :: u ! conserved variables
|
|
|
|
real, dimension(:,:,:,:,:), allocatable :: f ! numerical fluxes
|
2009-09-09 16:37:28 -03:00
|
|
|
end type block_data
|
|
|
|
|
2009-09-18 17:43:48 -03:00
|
|
|
! define block_info structure for boundary exchange
|
|
|
|
!
|
|
|
|
type block_info
|
|
|
|
type(block_info) , pointer :: prev ! pointer to the previous block
|
|
|
|
type(block_info) , pointer :: next ! pointer to the next block
|
|
|
|
type(block_meta) , pointer :: block ! pointer to the meta block
|
|
|
|
type(block_meta) , pointer :: neigh ! pointer to the neighbor block
|
2009-09-18 20:34:23 -03:00
|
|
|
integer(kind=4) :: direction ! direction of the neighbor block
|
|
|
|
integer(kind=4) :: side ! side of the neighbor block
|
|
|
|
integer(kind=4) :: face ! face of the neighbor block
|
|
|
|
integer(kind=4) :: level_difference ! the difference of levels
|
2009-09-18 17:43:48 -03:00
|
|
|
end type block_info
|
|
|
|
|
2009-09-09 16:37:28 -03:00
|
|
|
! chains of meta blocks and data blocks
|
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
type(block_meta), pointer, save :: list_meta, last_meta
|
|
|
|
type(block_data), pointer, save :: list_data, last_data
|
2009-09-09 16:37:28 -03:00
|
|
|
|
2008-11-11 16:12:26 -06:00
|
|
|
! stored last id (should always increase)
|
|
|
|
!
|
|
|
|
integer(kind=4) , save :: last_id
|
|
|
|
|
2009-05-18 22:46:19 +02:00
|
|
|
! store number of allocated blocks and leafs
|
2008-12-19 16:36:07 -06:00
|
|
|
!
|
2010-12-01 09:23:03 -02:00
|
|
|
integer(kind=4) , save :: mblocks, dblocks, nleafs
|
2008-12-19 16:36:07 -06:00
|
|
|
|
2011-05-05 15:56:38 -03:00
|
|
|
! store numbers of variables and fluxes
|
|
|
|
!
|
|
|
|
integer(kind=4) , save :: nvars, nflux
|
|
|
|
|
|
|
|
! store data block dimensions
|
|
|
|
!
|
|
|
|
integer(kind=4) , save :: nx, ny, nz
|
|
|
|
|
2008-11-04 21:00:50 -06:00
|
|
|
contains
|
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
! init_blocks: subroutine initializes the block variables
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
subroutine init_blocks()
|
2008-11-07 00:10:09 -06:00
|
|
|
|
2011-05-05 16:59:51 -03:00
|
|
|
use error, only : print_warning
|
2008-11-07 00:10:09 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2010-12-01 09:23:03 -02:00
|
|
|
! nullify list pointers
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
|
|
|
nullify(list_meta)
|
|
|
|
nullify(list_data)
|
2011-04-14 13:40:04 +02:00
|
|
|
nullify(last_meta)
|
|
|
|
nullify(last_data)
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2010-12-01 09:23:03 -02:00
|
|
|
! reset the number of meta and data blocks, and leafs
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2010-12-01 09:23:03 -02:00
|
|
|
mblocks = 0
|
2009-09-11 21:52:18 -03:00
|
|
|
dblocks = 0
|
2009-05-18 22:46:19 +02:00
|
|
|
nleafs = 0
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2011-05-05 15:56:38 -03:00
|
|
|
! set the initial number of variables and fluxes
|
|
|
|
!
|
|
|
|
nvars = 1
|
|
|
|
nflux = 0
|
|
|
|
|
|
|
|
! set the initial data block resolution
|
|
|
|
!
|
|
|
|
nx = 1
|
|
|
|
ny = 1
|
|
|
|
nz = 1
|
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
! reset ID
|
2008-12-05 14:24:01 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
last_id = 0
|
2011-04-15 00:25:22 +02:00
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
end subroutine init_blocks
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
! clear_blocks: subroutine clears the block variables
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2010-12-01 09:23:03 -02:00
|
|
|
subroutine clear_blocks()
|
2008-11-07 00:10:09 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
! pointers
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer :: pmeta
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2009-09-10 16:18:59 -03:00
|
|
|
! clear all meta blocks
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta => list_meta
|
|
|
|
do while(associated(pmeta))
|
|
|
|
call deallocate_metablock(pmeta)
|
|
|
|
pmeta => list_meta
|
2009-09-10 16:18:59 -03:00
|
|
|
end do
|
2010-10-11 22:46:07 -03:00
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2009-09-21 19:02:05 -03:00
|
|
|
end subroutine clear_blocks
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
! increase_id: function increases the last ID by 1 and returns it
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2010-09-18 11:42:09 +02:00
|
|
|
function increase_id()
|
2008-11-11 16:12:26 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
! return variable
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
integer(kind=4) :: increase_id
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-06 20:11:36 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
! increase ID by 1
|
2008-12-06 20:11:36 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
last_id = last_id + 1
|
2008-11-11 16:12:26 -06:00
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
! return ID
|
2008-12-06 20:11:36 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
increase_id = last_id
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
return
|
2008-11-04 21:00:50 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-05 22:16:24 -06:00
|
|
|
!
|
2009-09-21 19:02:05 -03:00
|
|
|
end function increase_id
|
2008-11-04 21:00:50 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-04 21:00:50 -06:00
|
|
|
!
|
2011-05-05 15:56:38 -03:00
|
|
|
! set_datablock_dims: subroutine set the number of variables and dimensions
|
|
|
|
! for arrays allocated in data blocks
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine set_datablock_dims(nv, nf, ni, nj, nk)
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
integer(kind=4), intent(in) :: nv, nf, ni, nj, nk
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
nvars = nv
|
|
|
|
nflux = nf
|
|
|
|
nx = ni
|
|
|
|
ny = nj
|
|
|
|
nz = nk
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine set_datablock_dims
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
! allocate_metablock: subroutine allocates space for one meta block and returns
|
|
|
|
! the pointer to this block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine allocate_metablock(pmeta)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(out) :: pmeta
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: i, j, k
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! allocate block structure
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
allocate(pmeta)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
! nullify pointers
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%prev)
|
|
|
|
nullify(pmeta%next)
|
|
|
|
nullify(pmeta%parent)
|
|
|
|
nullify(pmeta%data)
|
2009-09-10 17:25:28 -03:00
|
|
|
do i = 1, nchild
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%child(i)%ptr)
|
2009-09-10 17:25:28 -03:00
|
|
|
end do
|
2009-09-11 21:52:18 -03:00
|
|
|
do k = 1, nfaces
|
|
|
|
do j = 1, nsides
|
2009-09-10 17:25:28 -03:00
|
|
|
do i = 1, ndims
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%neigh(i,j,k)%ptr)
|
2009-09-10 17:25:28 -03:00
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
|
|
|
! set unique ID
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%id = increase_id()
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
! unset the CPU number of current block, level, the configuration, refine and
|
|
|
|
! leaf flags
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%cpu = -1
|
|
|
|
pmeta%level = -1
|
|
|
|
pmeta%config = -1
|
|
|
|
pmeta%refine = 0
|
|
|
|
pmeta%leaf = .false.
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2010-12-01 09:31:38 -02:00
|
|
|
! initialize the position in the parent block
|
|
|
|
!
|
|
|
|
pmeta%pos(:) = -1
|
|
|
|
|
2010-10-11 02:10:47 -03:00
|
|
|
! initialize the coordinate
|
|
|
|
!
|
|
|
|
pmeta%coord(:) = 0
|
|
|
|
|
2009-09-25 11:41:19 -03:00
|
|
|
! initialize bounds of the block
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%xmin = 0.0
|
|
|
|
pmeta%xmax = 1.0
|
|
|
|
pmeta%ymin = 0.0
|
|
|
|
pmeta%ymax = 1.0
|
|
|
|
pmeta%zmin = 0.0
|
|
|
|
pmeta%zmax = 1.0
|
2009-09-25 11:41:19 -03:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! increase the number of allocated meta blocks
|
|
|
|
!
|
2010-12-01 09:23:03 -02:00
|
|
|
mblocks = mblocks + 1
|
2009-09-21 19:02:05 -03:00
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine allocate_metablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! allocate_datablock: subroutine allocates space for one data block and returns
|
|
|
|
! the pointer to this block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine allocate_datablock(pdata)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_data), pointer, intent(out) :: pdata
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! allocate block structure
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
allocate(pdata)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
! nullify pointers
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pdata%prev)
|
|
|
|
nullify(pdata%next)
|
|
|
|
nullify(pdata%meta)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2010-07-26 20:46:45 -03:00
|
|
|
! allocate space for the conserved variables
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
2011-05-05 15:56:38 -03:00
|
|
|
allocate(pdata%u(nvars,nx,ny,nz))
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2011-04-30 12:11:57 -03:00
|
|
|
! allocate space for the numerical fluxes
|
|
|
|
!
|
2011-05-05 15:56:38 -03:00
|
|
|
if (nflux .gt. 0) allocate(pdata%f(NDIMS,nflux,nx,ny,nz))
|
2011-04-30 12:11:57 -03:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! increase the number of allocated meta blocks
|
|
|
|
!
|
|
|
|
dblocks = dblocks + 1
|
2008-11-04 21:00:50 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-05 22:16:24 -06:00
|
|
|
!
|
2009-09-21 19:02:05 -03:00
|
|
|
end subroutine allocate_datablock
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2009-09-10 16:18:59 -03:00
|
|
|
! deallocate_metablock: subroutine deallocates space ocuppied by a given metablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine deallocate_metablock(pmeta)
|
2009-09-10 16:18:59 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
2009-09-10 16:18:59 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: i, j, k
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(pmeta)) then
|
|
|
|
|
|
|
|
! if this is the first block in the list, update the list_meta pointer
|
|
|
|
!
|
|
|
|
if (pmeta%id .eq. list_meta%id) &
|
|
|
|
list_meta => pmeta%next
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2009-09-26 14:27:47 -03:00
|
|
|
! if this is the last block in the list, update the last_meta pointer
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (pmeta%id .eq. last_meta%id) &
|
|
|
|
last_meta => pmeta%prev
|
2009-09-10 16:18:59 -03:00
|
|
|
|
|
|
|
! update the pointer of previous and next blocks
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(pmeta%prev)) &
|
|
|
|
pmeta%prev%next => pmeta%next
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(pmeta%next)) &
|
|
|
|
pmeta%next%prev => pmeta%prev
|
2009-09-10 16:18:59 -03:00
|
|
|
|
|
|
|
! nullify children
|
|
|
|
!
|
|
|
|
do i = 1, nchild
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%child(i)%ptr)
|
2009-09-10 16:18:59 -03:00
|
|
|
end do
|
|
|
|
|
|
|
|
! nullify neighbors
|
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
do k = 1, nfaces
|
|
|
|
do j = 1, nsides
|
2009-09-10 16:18:59 -03:00
|
|
|
do i = 1, ndims
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%neigh(i,j,k)%ptr)
|
2009-09-10 16:18:59 -03:00
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
! if corresponding data block is allocated, deallocate it too
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(pmeta%data)) &
|
|
|
|
call deallocate_datablock(pmeta%data)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-09-10 16:18:59 -03:00
|
|
|
! nullify pointers
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%next)
|
|
|
|
nullify(pmeta%prev)
|
|
|
|
nullify(pmeta%data)
|
|
|
|
nullify(pmeta%parent)
|
2009-09-10 16:18:59 -03:00
|
|
|
|
|
|
|
! free and nullify the block
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
deallocate(pmeta)
|
|
|
|
nullify(pmeta)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
|
|
|
! decrease the number of allocated blocks
|
|
|
|
!
|
2010-12-01 09:23:03 -02:00
|
|
|
mblocks = mblocks - 1
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
|
|
|
!
|
2009-09-10 16:18:59 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine deallocate_metablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! deallocate_datablock: subroutine deallocates space ocuppied by a given datablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine deallocate_datablock(pdata)
|
2009-09-10 16:18:59 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_data), pointer, intent(inout) :: pdata
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(pdata)) then
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2009-09-25 19:24:02 -03:00
|
|
|
! if this is the first block in the list, update the list_data pointer
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (pdata%meta%id .eq. list_data%meta%id) &
|
|
|
|
list_data => pdata%next
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2009-09-25 19:24:02 -03:00
|
|
|
! if this is the last block in the list, update the last_data pointer
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (pdata%meta%id .eq. last_data%meta%id) &
|
|
|
|
last_data => pdata%prev
|
2009-09-25 19:24:02 -03:00
|
|
|
|
2009-09-10 16:18:59 -03:00
|
|
|
! update the pointer of previous and next blocks
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(pdata%prev)) &
|
|
|
|
pdata%prev%next => pdata%next
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(pdata%next)) &
|
|
|
|
pdata%next%prev => pdata%prev
|
2009-09-10 16:18:59 -03:00
|
|
|
|
|
|
|
! deallocate variables
|
|
|
|
!
|
2011-05-05 15:56:38 -03:00
|
|
|
if (allocated(pdata%u)) deallocate(pdata%u)
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2011-04-30 12:11:57 -03:00
|
|
|
! deallocate numerical fluxes
|
|
|
|
!
|
2011-05-05 15:56:38 -03:00
|
|
|
if (allocated(pdata%f)) deallocate(pdata%f)
|
2011-04-30 12:11:57 -03:00
|
|
|
|
2009-09-10 16:18:59 -03:00
|
|
|
! nullify pointers
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pdata%next)
|
|
|
|
nullify(pdata%prev)
|
|
|
|
nullify(pdata%meta)
|
2009-09-10 16:18:59 -03:00
|
|
|
|
|
|
|
! free and nullify the block
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
deallocate(pdata)
|
|
|
|
nullify(pdata)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
|
|
|
! decrease the number of allocated blocks
|
|
|
|
!
|
|
|
|
dblocks = dblocks - 1
|
|
|
|
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2009-09-21 19:02:05 -03:00
|
|
|
end subroutine deallocate_datablock
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
! append_metablock: subroutine allocates space for one meta block and appends it
|
|
|
|
! to the meta block list
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine append_metablock(pmeta)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(out) :: pmeta
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! allocate block
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
call allocate_metablock(pmeta)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
! add to the list
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(last_meta)) then
|
|
|
|
pmeta%prev => last_meta
|
|
|
|
last_meta%next => pmeta
|
2011-04-14 13:40:04 +02:00
|
|
|
else
|
|
|
|
list_meta => pmeta
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
! set the pointer to the last block in the list
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
last_meta => pmeta
|
2009-09-21 19:02:05 -03:00
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine append_metablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! append_datablock: subroutine allocates space for one data block and appends it
|
|
|
|
! to the data block list
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine append_datablock(pdata)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_data), pointer, intent(out) :: pdata
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! allocate block
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
call allocate_datablock(pdata)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
! add to the list
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(last_data)) then
|
|
|
|
pdata%prev => last_data
|
|
|
|
last_data%next => pdata
|
2011-04-14 13:40:04 +02:00
|
|
|
else
|
|
|
|
list_data => pdata
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
! set the pointer to the last block in the list
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
last_data => pdata
|
2009-09-21 19:02:05 -03:00
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine append_datablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-10 17:46:36 -03:00
|
|
|
! associate_blocks: subroutine associates a pair of meta and data blocks
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine associate_blocks(pmeta, pdata)
|
2009-09-10 17:46:36 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
type(block_data), pointer, intent(inout) :: pdata
|
2009-09-10 17:46:36 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%data => pdata
|
|
|
|
pdata%meta => pmeta
|
2009-09-10 17:46:36 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine associate_blocks
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-10 18:15:30 -03:00
|
|
|
! metablock_setleaf: subroutine sets the leaf flag of data block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine metablock_setleaf(pmeta)
|
2009-09-10 18:15:30 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input/output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set the leaf flag
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%leaf = .true.
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
! increase the number of leafs
|
|
|
|
!
|
|
|
|
nleafs = nleafs + 1
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_setleaf
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! metablock_unsetleaf: subroutine unsets the leaf flag of data block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine metablock_unsetleaf(pmeta)
|
2009-09-10 18:15:30 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input/output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set the leaf flag
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%leaf = .false.
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
! decrease the number of leafs
|
|
|
|
!
|
|
|
|
nleafs = nleafs - 1
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_unsetleaf
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! metablock_setconfig: subroutine sets the config flag of data block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine metablock_setconfig(pmeta, config)
|
2009-09-10 18:15:30 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input/output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
2009-09-10 18:15:30 -03:00
|
|
|
integer(kind=4) , intent(in) :: config
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set the config flag
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%config = config
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_setconfig
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! metablock_setlevel: subroutine sets the level of data block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine metablock_setlevel(pmeta, level)
|
2009-09-10 18:15:30 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input/output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
2009-09-10 18:15:30 -03:00
|
|
|
integer(kind=4) , intent(in) :: level
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set the refinement level
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%level = level
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_setlevel
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2010-12-01 09:31:38 -02:00
|
|
|
! metablock_set_position: subroutine sets the position of the meta block in the
|
|
|
|
! parent block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine metablock_set_position(pmeta, px, py, pz)
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input/output arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
integer(kind=4) , intent(in) :: px, py, pz
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set the position in the parent block
|
|
|
|
!
|
|
|
|
pmeta%pos(1) = px
|
|
|
|
pmeta%pos(2) = py
|
|
|
|
#if NDIMS == 3
|
|
|
|
pmeta%pos(3) = pz
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_set_position
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2010-10-11 02:10:47 -03:00
|
|
|
! metablock_set_coord: subroutine sets the coordinates of the meta block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine metablock_set_coord(pmeta, px, py, pz)
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input/output arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
integer(kind=4) , intent(in) :: px, py, pz
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set the coordintaes
|
|
|
|
!
|
|
|
|
pmeta%coord(1) = px
|
|
|
|
pmeta%coord(2) = py
|
|
|
|
#if NDIMS == 3
|
|
|
|
pmeta%coord(3) = pz
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_set_coord
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-25 11:41:19 -03:00
|
|
|
! metablock_setbounds: subroutine sets the bounds of data block
|
2009-09-10 17:46:36 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine metablock_setbounds(pmeta, xmin, xmax, ymin, ymax, zmin, zmax)
|
2009-09-10 17:46:36 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input/output arguments
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
2009-09-10 18:15:30 -03:00
|
|
|
real , intent(in) :: xmin, xmax, ymin, ymax, zmin, zmax
|
2009-09-10 17:46:36 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set bounds of the block
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
pmeta%xmin = xmin
|
|
|
|
pmeta%xmax = xmax
|
|
|
|
pmeta%ymin = ymin
|
|
|
|
pmeta%ymax = ymax
|
|
|
|
pmeta%zmin = zmin
|
|
|
|
pmeta%zmax = zmax
|
2009-09-10 17:46:36 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2009-09-25 11:41:19 -03:00
|
|
|
end subroutine metablock_setbounds
|
2009-09-10 17:46:36 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2008-12-05 14:04:12 -06:00
|
|
|
! refine_block: subroutine refines selected block
|
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2011-05-05 16:51:35 -03:00
|
|
|
subroutine refine_block(pblock, res, falloc_data)
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2011-05-05 15:56:38 -03:00
|
|
|
use error, only : print_error
|
2008-12-05 14:04:12 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input parameters
|
|
|
|
!
|
2011-05-05 16:51:35 -03:00
|
|
|
type(block_meta), pointer , intent(inout) :: pblock
|
|
|
|
integer(kind=4), dimension(3), intent(in) :: res
|
|
|
|
logical , intent(in) :: falloc_data
|
2008-12-05 14:04:12 -06:00
|
|
|
|
|
|
|
! pointers
|
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
type(block_meta), pointer :: pneigh, pchild, pfirst, plast
|
|
|
|
type(block_data), pointer :: pdata
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
2009-09-13 22:58:55 -03:00
|
|
|
integer, dimension(nchild) :: config, order
|
|
|
|
integer, dimension(ndims,nsides,nfaces) :: set
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2010-10-11 22:46:07 -03:00
|
|
|
integer :: p, i, j, k, ic, jc, kc
|
2009-09-11 21:52:18 -03:00
|
|
|
real :: xln, yln, zln, xmn, xmx, ymn, ymx, zmn, zmx
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
|
|
|
! check if pointer is associated
|
|
|
|
!
|
|
|
|
if (associated(pblock)) then
|
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! unset block leaf flag
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
call metablock_unsetleaf(pblock)
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! reset refinement flag
|
2009-05-18 22:46:19 +02:00
|
|
|
!
|
|
|
|
pblock%refine = 0
|
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! iterate over all child blocks
|
2009-05-18 22:46:19 +02:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
do p = 1, nchild
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! create child meta and data blocks
|
2009-05-18 22:46:19 +02:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
call allocate_metablock(pblock%child(p)%ptr)
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! set it as a leaf
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
call metablock_setleaf(pblock%child(p)%ptr)
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! assign pointer to the parent block
|
2008-12-05 14:22:02 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pblock%child(p)%ptr%parent => pblock
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! increase the refinement level
|
2009-01-02 20:18:57 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pblock%child(p)%ptr%level = pblock%level + 1
|
2009-01-02 20:18:57 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! copy the parent cpu number to each child
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pblock%child(p)%ptr%cpu = pblock%cpu
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
end do
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! assign neighbors of the child blocks
|
2009-09-13 22:58:55 -03:00
|
|
|
!
|
|
|
|
! interior of the block
|
2009-09-11 21:52:18 -03:00
|
|
|
!
|
|
|
|
do p = 1, nfaces
|
|
|
|
|
|
|
|
! X direction (left side)
|
|
|
|
!
|
|
|
|
pblock%child(2)%ptr%neigh(1,1,p)%ptr => pblock%child(1)%ptr
|
|
|
|
pblock%child(4)%ptr%neigh(1,1,p)%ptr => pblock%child(3)%ptr
|
|
|
|
#if NDIMS == 3
|
|
|
|
pblock%child(6)%ptr%neigh(1,1,p)%ptr => pblock%child(5)%ptr
|
|
|
|
pblock%child(8)%ptr%neigh(1,1,p)%ptr => pblock%child(7)%ptr
|
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-13 22:58:55 -03:00
|
|
|
pneigh => pblock%neigh(1,1,1)%ptr
|
|
|
|
if (associated(pneigh)) then
|
|
|
|
if (pneigh%id .eq. pblock%id) then
|
|
|
|
pblock%child(1)%ptr%neigh(1,1,p)%ptr => pblock%child(2)%ptr
|
|
|
|
pblock%child(3)%ptr%neigh(1,1,p)%ptr => pblock%child(4)%ptr
|
|
|
|
#if NDIMS == 3
|
|
|
|
pblock%child(5)%ptr%neigh(1,1,p)%ptr => pblock%child(6)%ptr
|
|
|
|
pblock%child(7)%ptr%neigh(1,1,p)%ptr => pblock%child(8)%ptr
|
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
|
|
|
end if
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
! X direction (right side)
|
|
|
|
!
|
|
|
|
pblock%child(1)%ptr%neigh(1,2,p)%ptr => pblock%child(2)%ptr
|
|
|
|
pblock%child(3)%ptr%neigh(1,2,p)%ptr => pblock%child(4)%ptr
|
|
|
|
#if NDIMS == 3
|
|
|
|
pblock%child(5)%ptr%neigh(1,2,p)%ptr => pblock%child(6)%ptr
|
|
|
|
pblock%child(7)%ptr%neigh(1,2,p)%ptr => pblock%child(8)%ptr
|
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-13 22:58:55 -03:00
|
|
|
pneigh => pblock%neigh(1,2,1)%ptr
|
|
|
|
if (associated(pneigh)) then
|
|
|
|
if (pneigh%id .eq. pblock%id) then
|
|
|
|
pblock%child(2)%ptr%neigh(1,2,p)%ptr => pblock%child(1)%ptr
|
|
|
|
pblock%child(4)%ptr%neigh(1,2,p)%ptr => pblock%child(3)%ptr
|
|
|
|
#if NDIMS == 3
|
|
|
|
pblock%child(6)%ptr%neigh(1,2,p)%ptr => pblock%child(5)%ptr
|
|
|
|
pblock%child(8)%ptr%neigh(1,2,p)%ptr => pblock%child(7)%ptr
|
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
|
|
|
end if
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
! Y direction (left side)
|
|
|
|
!
|
|
|
|
pblock%child(3)%ptr%neigh(2,1,p)%ptr => pblock%child(1)%ptr
|
|
|
|
pblock%child(4)%ptr%neigh(2,1,p)%ptr => pblock%child(2)%ptr
|
|
|
|
#if NDIMS == 3
|
|
|
|
pblock%child(7)%ptr%neigh(2,1,p)%ptr => pblock%child(5)%ptr
|
|
|
|
pblock%child(8)%ptr%neigh(2,1,p)%ptr => pblock%child(6)%ptr
|
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-13 22:58:55 -03:00
|
|
|
pneigh => pblock%neigh(2,1,1)%ptr
|
|
|
|
if (associated(pneigh)) then
|
|
|
|
if (pneigh%id .eq. pblock%id) then
|
|
|
|
pblock%child(1)%ptr%neigh(2,1,p)%ptr => pblock%child(3)%ptr
|
|
|
|
pblock%child(2)%ptr%neigh(2,1,p)%ptr => pblock%child(4)%ptr
|
|
|
|
#if NDIMS == 3
|
|
|
|
pblock%child(5)%ptr%neigh(2,1,p)%ptr => pblock%child(7)%ptr
|
|
|
|
pblock%child(6)%ptr%neigh(2,1,p)%ptr => pblock%child(8)%ptr
|
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
|
|
|
end if
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
! Y direction (right side)
|
|
|
|
!
|
|
|
|
pblock%child(1)%ptr%neigh(2,2,p)%ptr => pblock%child(3)%ptr
|
|
|
|
pblock%child(2)%ptr%neigh(2,2,p)%ptr => pblock%child(4)%ptr
|
|
|
|
#if NDIMS == 3
|
|
|
|
pblock%child(5)%ptr%neigh(2,2,p)%ptr => pblock%child(7)%ptr
|
|
|
|
pblock%child(6)%ptr%neigh(2,2,p)%ptr => pblock%child(8)%ptr
|
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-13 22:58:55 -03:00
|
|
|
pneigh => pblock%neigh(2,2,1)%ptr
|
|
|
|
if (associated(pneigh)) then
|
|
|
|
if (pneigh%id .eq. pblock%id) then
|
|
|
|
pblock%child(3)%ptr%neigh(2,2,p)%ptr => pblock%child(1)%ptr
|
|
|
|
pblock%child(4)%ptr%neigh(2,2,p)%ptr => pblock%child(2)%ptr
|
|
|
|
#if NDIMS == 3
|
|
|
|
pblock%child(7)%ptr%neigh(2,2,p)%ptr => pblock%child(5)%ptr
|
|
|
|
pblock%child(8)%ptr%neigh(2,2,p)%ptr => pblock%child(6)%ptr
|
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
|
|
|
end if
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
#if NDIMS == 3
|
|
|
|
! Z direction (left side)
|
|
|
|
!
|
|
|
|
pblock%child(5)%ptr%neigh(3,1,p)%ptr => pblock%child(1)%ptr
|
|
|
|
pblock%child(6)%ptr%neigh(3,1,p)%ptr => pblock%child(2)%ptr
|
|
|
|
pblock%child(7)%ptr%neigh(3,1,p)%ptr => pblock%child(3)%ptr
|
|
|
|
pblock%child(8)%ptr%neigh(3,1,p)%ptr => pblock%child(4)%ptr
|
2009-09-13 22:58:55 -03:00
|
|
|
pneigh => pblock%neigh(3,1,1)%ptr
|
|
|
|
if (associated(pneigh)) then
|
|
|
|
if (pneigh%id .eq. pblock%id) then
|
|
|
|
pblock%child(1)%ptr%neigh(3,1,p)%ptr => pblock%child(5)%ptr
|
|
|
|
pblock%child(2)%ptr%neigh(3,1,p)%ptr => pblock%child(6)%ptr
|
|
|
|
pblock%child(3)%ptr%neigh(3,1,p)%ptr => pblock%child(7)%ptr
|
|
|
|
pblock%child(4)%ptr%neigh(3,1,p)%ptr => pblock%child(8)%ptr
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
|
|
|
end if
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
! Z direction (right side)
|
|
|
|
!
|
|
|
|
pblock%child(1)%ptr%neigh(3,2,p)%ptr => pblock%child(5)%ptr
|
|
|
|
pblock%child(2)%ptr%neigh(3,2,p)%ptr => pblock%child(6)%ptr
|
|
|
|
pblock%child(3)%ptr%neigh(3,2,p)%ptr => pblock%child(7)%ptr
|
|
|
|
pblock%child(4)%ptr%neigh(3,2,p)%ptr => pblock%child(8)%ptr
|
2009-09-13 22:58:55 -03:00
|
|
|
pneigh => pblock%neigh(3,2,1)%ptr
|
|
|
|
if (associated(pneigh)) then
|
|
|
|
if (pneigh%id .eq. pblock%id) then
|
|
|
|
pblock%child(5)%ptr%neigh(3,2,p)%ptr => pblock%child(1)%ptr
|
|
|
|
pblock%child(6)%ptr%neigh(3,2,p)%ptr => pblock%child(2)%ptr
|
|
|
|
pblock%child(7)%ptr%neigh(3,2,p)%ptr => pblock%child(3)%ptr
|
|
|
|
pblock%child(8)%ptr%neigh(3,2,p)%ptr => pblock%child(4)%ptr
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
|
|
|
end if
|
2009-09-11 21:52:18 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
end do
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
! prepare set array
|
2008-12-06 20:11:36 -06:00
|
|
|
!
|
2009-09-13 22:58:55 -03:00
|
|
|
#if NDIMS == 2
|
|
|
|
set(1,1,:) = (/ 1, 3 /)
|
|
|
|
set(1,2,:) = (/ 2, 4 /)
|
|
|
|
set(2,1,:) = (/ 1, 2 /)
|
|
|
|
set(2,2,:) = (/ 3, 4 /)
|
|
|
|
#endif /* NDIMS == 2 */
|
2009-09-11 21:52:18 -03:00
|
|
|
#if NDIMS == 3
|
2009-09-13 22:58:55 -03:00
|
|
|
set(1,1,:) = (/ 1, 3, 5, 7 /)
|
|
|
|
set(1,2,:) = (/ 2, 4, 6, 8 /)
|
|
|
|
set(2,1,:) = (/ 1, 2, 5, 6 /)
|
|
|
|
set(2,2,:) = (/ 3, 4, 7, 8 /)
|
|
|
|
set(3,1,:) = (/ 1, 2, 3, 4 /)
|
|
|
|
set(3,2,:) = (/ 5, 6, 7, 8 /)
|
2009-09-11 21:52:18 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
! set pointers to neighbors and update neighbors pointers
|
2009-09-11 21:52:18 -03:00
|
|
|
!
|
2009-09-13 22:58:55 -03:00
|
|
|
do i = 1, ndims
|
|
|
|
do j = 1, nsides
|
|
|
|
do k = 1, nfaces
|
|
|
|
pneigh => pblock%neigh(i,j,k)%ptr
|
2011-04-23 17:17:11 -03:00
|
|
|
pchild => pblock%child(set(i,j,k))%ptr
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
if (associated(pneigh)) then
|
|
|
|
if (pneigh%id .ne. pblock%id) then
|
2008-12-19 16:36:07 -06:00
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
! point to the right neighbor
|
2009-09-11 21:52:18 -03:00
|
|
|
!
|
2009-09-13 22:58:55 -03:00
|
|
|
do p = 1, nfaces
|
2011-04-23 17:17:11 -03:00
|
|
|
pchild%neigh(i,j,p)%ptr => pneigh
|
2009-09-13 22:58:55 -03:00
|
|
|
end do
|
2008-12-19 16:36:07 -06:00
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
! neighbor level is the same as the refined block
|
2009-09-11 21:52:18 -03:00
|
|
|
!
|
2009-09-13 22:58:55 -03:00
|
|
|
if (pneigh%level .eq. pblock%level) then
|
2011-04-23 17:17:11 -03:00
|
|
|
pneigh%neigh(i,3-j,k)%ptr => pchild
|
2009-09-13 22:58:55 -03:00
|
|
|
end if
|
2008-12-19 16:36:07 -06:00
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
! neighbor level is the same as the child block
|
2009-09-11 21:52:18 -03:00
|
|
|
!
|
2009-09-13 22:58:55 -03:00
|
|
|
if (pneigh%level .gt. pblock%level) then
|
|
|
|
do p = 1, nfaces
|
2011-04-23 17:17:11 -03:00
|
|
|
pneigh%neigh(i,3-j,p)%ptr => pchild
|
2009-09-13 22:58:55 -03:00
|
|
|
end do
|
|
|
|
end if
|
2008-12-19 16:36:07 -06:00
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
end if
|
2008-12-19 16:36:07 -06:00
|
|
|
|
2009-09-13 22:58:55 -03:00
|
|
|
end if
|
2008-12-19 16:36:07 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
end do
|
2009-09-13 22:58:55 -03:00
|
|
|
end do
|
|
|
|
end do
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2011-04-23 21:48:51 -03:00
|
|
|
! reset neighbor pointers of the parent block
|
|
|
|
!
|
|
|
|
do i = 1, ndims
|
|
|
|
do j = 1, nsides
|
|
|
|
do k = 1, nfaces
|
|
|
|
nullify(pblock%neigh(i,j,k)%ptr)
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! set corresponding configuration of the new blocks
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
|
|
|
select case(pblock%config)
|
2009-10-03 20:11:03 -03:00
|
|
|
case(0)
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
#if NDIMS == 2
|
|
|
|
config(:) = (/ 0, 0, 0, 0 /)
|
|
|
|
order (:) = (/ 1, 2, 3, 4 /)
|
|
|
|
#endif /* NDIMS == 2 */
|
|
|
|
#if NDIMS == 3
|
|
|
|
config(:) = (/ 0, 0, 0, 0, 0, 0, 0, 0 /)
|
|
|
|
order (:) = (/ 1, 2, 3, 4, 5, 6, 7, 8 /)
|
|
|
|
#endif /* NDIMS == 3 */
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
case(12)
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
#if NDIMS == 2
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 13, 12, 12, 42 /)
|
|
|
|
order (:) = (/ 1, 3, 4, 2 /)
|
2009-09-11 21:52:18 -03:00
|
|
|
#endif /* NDIMS == 2 */
|
2009-09-13 22:58:55 -03:00
|
|
|
#if NDIMS == 3
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 13, 15, 15, 78, 78, 62, 62, 42 /)
|
|
|
|
order (:) = (/ 1, 3, 7, 5, 6, 8, 4, 2 /)
|
2009-09-13 22:58:55 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
case(13)
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
#if NDIMS == 2
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 12, 13, 13, 43 /)
|
|
|
|
order (:) = (/ 1, 2, 4, 3 /)
|
2009-09-11 21:52:18 -03:00
|
|
|
#endif /* NDIMS == 2 */
|
2009-09-13 22:58:55 -03:00
|
|
|
#if NDIMS == 3
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 15, 12, 12, 68, 68, 43, 43, 73 /)
|
|
|
|
order (:) = (/ 1, 5, 6, 2, 4, 8, 7, 3 /)
|
|
|
|
|
|
|
|
case(15)
|
|
|
|
|
|
|
|
config(:) = (/ 12, 13, 13, 48, 48, 75, 75, 65 /)
|
|
|
|
order (:) = (/ 1, 2, 4, 3, 7, 8, 6, 5 /)
|
2009-09-13 22:58:55 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
case(42)
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
#if NDIMS == 2
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 43, 42, 42, 12 /)
|
|
|
|
order (:) = (/ 4, 3, 1, 2 /)
|
2009-09-11 21:52:18 -03:00
|
|
|
#endif /* NDIMS == 2 */
|
2009-09-13 22:58:55 -03:00
|
|
|
#if NDIMS == 3
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 48, 43, 43, 75, 75, 12, 12, 62 /)
|
|
|
|
order (:) = (/ 4, 8, 7, 3, 1, 5, 6, 2 /)
|
2009-09-13 22:58:55 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-11 21:52:18 -03:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
case(43)
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
#if NDIMS == 2
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 42, 43, 43, 13 /)
|
|
|
|
order (:) = (/ 4, 2, 1, 3 /)
|
2009-09-11 21:52:18 -03:00
|
|
|
#endif /* NDIMS == 2 */
|
2009-09-13 22:58:55 -03:00
|
|
|
#if NDIMS == 3
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 42, 48, 48, 65, 65, 73, 73, 13 /)
|
|
|
|
order (:) = (/ 4, 2, 6, 8, 7, 5, 1, 3 /)
|
2009-09-13 22:58:55 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
#if NDIMS == 3
|
2009-10-03 20:11:03 -03:00
|
|
|
case(48)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 43, 42, 42, 15, 15, 68, 68, 78 /)
|
|
|
|
order (:) = (/ 4, 3, 1, 2, 6, 5, 7, 8 /)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
case(62)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 65, 68, 68, 73, 73, 42, 42, 12 /)
|
|
|
|
order (:) = (/ 6, 5, 7, 8, 4, 3, 1, 2 /)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
case(65)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 68, 62, 62, 43, 43, 15, 15, 75 /)
|
|
|
|
order (:) = (/ 6, 8, 4, 2, 1, 3, 7, 5 /)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
case(68)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2009-10-03 20:11:03 -03:00
|
|
|
config(:) = (/ 62, 65, 65, 13, 13, 78, 78, 48 /)
|
|
|
|
order (:) = (/ 6, 2, 1, 5, 7, 3, 4 , 8 /)
|
|
|
|
|
|
|
|
case(73)
|
|
|
|
|
|
|
|
config(:) = (/ 78, 75, 75, 62, 62, 13, 13, 43 /)
|
|
|
|
order (:) = (/ 7, 8, 6, 5, 1, 2, 4, 3 /)
|
|
|
|
|
|
|
|
case(75)
|
|
|
|
|
|
|
|
config(:) = (/ 73, 78, 78, 42, 42, 65, 65, 15 /)
|
|
|
|
order (:) = (/ 7, 3, 4, 8, 6, 2, 1, 5 /)
|
|
|
|
|
|
|
|
case(78)
|
|
|
|
|
|
|
|
config(:) = (/ 75, 73, 73, 12, 12, 48, 48, 68 /)
|
|
|
|
order (:) = (/ 7, 5, 1, 3, 4, 2, 6, 8 /)
|
2009-09-13 22:58:55 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2009-09-11 21:52:18 -03:00
|
|
|
|
|
|
|
end select
|
2008-12-05 14:04:12 -06:00
|
|
|
|
|
|
|
! set blocks configurations
|
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
do p = 1, nchild
|
2009-09-13 22:58:55 -03:00
|
|
|
pblock%child(order(p))%ptr%config = config(p)
|
2009-09-11 21:52:18 -03:00
|
|
|
end do
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! connect blocks in chain
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
do p = 2, nchild
|
|
|
|
pblock%child(order(p ))%ptr%prev => pblock%child(order(p-1))%ptr
|
|
|
|
pblock%child(order(p-1))%ptr%next => pblock%child(order(p ))%ptr
|
|
|
|
end do
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! insert this chain after the parent block
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pneigh => pblock%next
|
|
|
|
pfirst => pblock%child(order( 1))%ptr
|
|
|
|
plast => pblock%child(order(nchild))%ptr
|
|
|
|
if (associated(pneigh)) then
|
|
|
|
pneigh%prev => plast
|
|
|
|
plast%next => pneigh
|
|
|
|
else
|
|
|
|
last_meta => plast
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
pblock%next => pfirst
|
|
|
|
pfirst%prev => pblock
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! calculate the size of new blocks
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-25 11:41:19 -03:00
|
|
|
xln = 0.5 * (pblock%xmax - pblock%xmin)
|
|
|
|
yln = 0.5 * (pblock%ymax - pblock%ymin)
|
2009-09-11 21:52:18 -03:00
|
|
|
#if NDIMS == 3
|
2009-09-25 11:41:19 -03:00
|
|
|
zln = 0.5 * (pblock%zmax - pblock%zmin)
|
2009-09-11 21:52:18 -03:00
|
|
|
#else /* NDIMS == 3 */
|
2009-09-25 11:41:19 -03:00
|
|
|
zln = (pblock%zmax - pblock%zmin)
|
2009-09-11 21:52:18 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! iterate over all children and allocate data blocks
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-25 11:41:19 -03:00
|
|
|
do p = 1, nchild
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! assign a pointer to the current child
|
|
|
|
!
|
2009-09-25 11:41:19 -03:00
|
|
|
pchild => pblock%child(p)%ptr
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2010-10-11 22:46:07 -03:00
|
|
|
! calculate the block position indices
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-25 11:41:19 -03:00
|
|
|
i = mod((p - 1) ,2)
|
|
|
|
j = mod((p - 1) / 2,2)
|
|
|
|
k = mod((p - 1) / 4,2)
|
2009-09-11 21:52:18 -03:00
|
|
|
|
2010-12-01 09:31:38 -02:00
|
|
|
! set the block position
|
|
|
|
!
|
|
|
|
call metablock_set_position(pchild, i, j, k)
|
|
|
|
|
2010-10-11 22:46:07 -03:00
|
|
|
! set the block coordinates
|
|
|
|
!
|
2011-05-05 16:51:35 -03:00
|
|
|
ic = pblock%coord(1) + i * res(1)
|
|
|
|
jc = pblock%coord(2) + j * res(2)
|
2010-10-11 22:46:07 -03:00
|
|
|
#if NDIMS == 3
|
2011-05-05 16:51:35 -03:00
|
|
|
kc = pblock%coord(3) + k * res(3)
|
2010-10-11 22:46:07 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
call metablock_set_coord(pchild, ic, jc, kc)
|
|
|
|
|
|
|
|
! calculate block bounds
|
|
|
|
!
|
2009-09-25 11:41:19 -03:00
|
|
|
xmn = pblock%xmin + xln * i
|
|
|
|
ymn = pblock%ymin + yln * j
|
|
|
|
zmn = pblock%zmin + zln * k
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-25 11:41:19 -03:00
|
|
|
xmx = xmn + xln
|
|
|
|
ymx = ymn + yln
|
|
|
|
zmx = zmn + zln
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! set block bounds
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-25 11:41:19 -03:00
|
|
|
call metablock_setbounds(pchild, xmn, xmx, ymn, ymx, zmn, zmx)
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! allocate data blocks if necessary
|
|
|
|
!
|
|
|
|
if (falloc_data) then
|
|
|
|
|
|
|
|
! iterate over all children and allocate data blocks
|
|
|
|
!
|
|
|
|
do p = 1, nchild
|
|
|
|
|
|
|
|
! assign a pointer to the current child
|
|
|
|
!
|
|
|
|
pchild => pblock%child(p)%ptr
|
|
|
|
|
|
|
|
! allocate data block
|
|
|
|
!
|
|
|
|
call allocate_datablock(pdata)
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! associate with the meta block
|
|
|
|
!
|
|
|
|
call associate_blocks(pchild, pdata)
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
end do
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! connect blocks in chain
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
do p = 2, nchild
|
|
|
|
pblock%child(order(p ))%ptr%data%prev => pblock%child(order(p-1))%ptr%data
|
|
|
|
pblock%child(order(p-1))%ptr%data%next => pblock%child(order(p ))%ptr%data
|
|
|
|
end do
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! insert this chain after the parent block
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pdata => pblock%data%next
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
pfirst => pblock%child(order( 1))%ptr
|
|
|
|
plast => pblock%child(order(nchild))%ptr
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
if (associated(pdata)) then
|
|
|
|
pdata%prev => plast%data
|
|
|
|
plast%data%next => pdata
|
|
|
|
else
|
|
|
|
last_data => plast%data
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
2008-12-05 14:04:12 -06:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
pblock%data%next => pfirst%data
|
|
|
|
pfirst%data%prev => pblock%data
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
! point the current block to the last created one
|
|
|
|
!
|
|
|
|
pblock => plast
|
|
|
|
|
2008-12-05 14:04:12 -06:00
|
|
|
else
|
|
|
|
|
|
|
|
! terminate program if the pointer passed by argument is not associated
|
|
|
|
!
|
2011-05-05 16:51:35 -03:00
|
|
|
call print_error("blocks::refine_block","Input pointer is not associated! Terminating!")
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
|
|
|
end subroutine refine_block
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
|
|
|
! derefine_block: subroutine derefines selected block
|
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
|
|
|
subroutine derefine_block(pblock)
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input parameters
|
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pblock
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
integer :: i, j, k, l, p
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-29 17:19:26 -03:00
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
integer, dimension(ndims, nsides, nfaces) :: arr
|
|
|
|
|
|
|
|
! local pointers
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
type(block_meta), pointer :: pchild, pneigh
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-29 17:19:26 -03:00
|
|
|
! prepare reference array
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
#if NDIMS == 3
|
2009-09-29 17:19:26 -03:00
|
|
|
arr(1,1,:) = (/ 1, 3, 5, 7 /)
|
|
|
|
arr(1,2,:) = (/ 2, 4, 6, 8 /)
|
|
|
|
arr(2,1,:) = (/ 1, 2, 5, 6 /)
|
|
|
|
arr(2,2,:) = (/ 3, 4, 7, 8 /)
|
|
|
|
arr(3,1,:) = (/ 1, 2, 3, 4 /)
|
|
|
|
arr(3,2,:) = (/ 5, 6, 7, 8 /)
|
|
|
|
#else /* NDIMS == 3 */
|
|
|
|
arr(1,1,:) = (/ 1, 3 /)
|
|
|
|
arr(1,2,:) = (/ 2, 4 /)
|
|
|
|
arr(2,1,:) = (/ 1, 2 /)
|
|
|
|
arr(2,2,:) = (/ 3, 4 /)
|
2009-09-21 16:57:34 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-29 17:19:26 -03:00
|
|
|
! iterate over all boundaries of the parent block
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
do i = 1, ndims
|
|
|
|
do j = 1, nsides
|
|
|
|
do k = 1, nfaces
|
2009-09-29 17:19:26 -03:00
|
|
|
|
|
|
|
! calculate the right child number
|
|
|
|
!
|
|
|
|
p = arr(i,j,k)
|
|
|
|
|
|
|
|
! assign the pointer to the current neighbor
|
|
|
|
!
|
|
|
|
pneigh => pblock%child(p)%ptr%neigh(i,j,k)%ptr
|
|
|
|
|
|
|
|
! assign the right neighbor to the current neighbor pointer
|
|
|
|
!
|
|
|
|
pblock%neigh(i,j,k)%ptr => pneigh
|
|
|
|
|
|
|
|
! update the neighbor fields of neighbors
|
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
if (associated(pneigh)) then
|
|
|
|
l = 3 - j
|
|
|
|
do p = 1, nfaces
|
|
|
|
pneigh%neigh(i,l,p)%ptr => pblock
|
|
|
|
end do
|
|
|
|
end if
|
2009-09-29 17:19:26 -03:00
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
2008-12-19 16:36:07 -06:00
|
|
|
|
2009-09-29 17:19:26 -03:00
|
|
|
! deallocate child blocks
|
|
|
|
!
|
|
|
|
do p = 1, nchild
|
|
|
|
call metablock_unsetleaf(pblock%child(p)%ptr)
|
|
|
|
call deallocate_metablock(pblock%child(p)%ptr)
|
|
|
|
end do
|
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
! set the leaf flag of parent block
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-21 17:16:29 -03:00
|
|
|
call metablock_setleaf(pblock)
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
! reset the refinement flag of the parent block
|
|
|
|
!
|
2008-12-13 15:08:18 -06:00
|
|
|
pblock%refine = 0
|
2009-09-21 16:57:34 -03:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
|
|
|
end subroutine derefine_block
|
2011-04-21 14:13:25 -03:00
|
|
|
#ifdef DEBUG
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! check_metablock: subroutine checks if the meta block has proper structure
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine check_metablock(pblock, string)
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input parameters
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(in) :: pblock
|
|
|
|
character(len=*) , intent(in) :: string
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: p, i, j, k
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_meta), pointer :: ptemp
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! check block ID
|
|
|
|
!
|
|
|
|
ptemp => pblock
|
|
|
|
if (ptemp%id .le. 0 .or. ptemp%id .gt. last_id) then
|
|
|
|
print *, ''
|
|
|
|
print *, ''
|
|
|
|
print *, trim(string)
|
|
|
|
print *, 'wrong meta block id = ', ptemp%id
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
|
|
|
|
! check prev ID
|
|
|
|
!
|
|
|
|
ptemp => pblock%prev
|
|
|
|
if (associated(ptemp)) then
|
|
|
|
if (ptemp%id .le. 0 .or. ptemp%id .gt. last_id) then
|
|
|
|
print *, ''
|
|
|
|
print *, ''
|
|
|
|
print *, trim(string)
|
|
|
|
print *, 'wrong previous block id = ', ptemp%id, pblock%id
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
! check next ID
|
|
|
|
!
|
|
|
|
ptemp => pblock%next
|
|
|
|
if (associated(ptemp)) then
|
|
|
|
if (ptemp%id .le. 0 .or. ptemp%id .gt. last_id) then
|
|
|
|
print *, ''
|
|
|
|
print *, ''
|
|
|
|
print *, trim(string)
|
|
|
|
print *, 'wrong next block id = ', ptemp%id, pblock%id
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
! check parent ID
|
|
|
|
!
|
|
|
|
ptemp => pblock%parent
|
|
|
|
if (associated(ptemp)) then
|
|
|
|
if (ptemp%id .le. 0 .or. ptemp%id .gt. last_id) then
|
|
|
|
print *, ''
|
|
|
|
print *, ''
|
|
|
|
print *, trim(string)
|
|
|
|
print *, 'wrong parent block id = ', ptemp%id, pblock%id
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
! check children IDs
|
|
|
|
!
|
|
|
|
do p = 1, nchild
|
|
|
|
ptemp => pblock%child(p)%ptr
|
|
|
|
if (associated(ptemp)) then
|
|
|
|
if (ptemp%id .le. 0 .or. ptemp%id .gt. last_id) then
|
|
|
|
print *, ''
|
|
|
|
print *, ''
|
|
|
|
print *, trim(string)
|
|
|
|
print *, 'wrong child block id = ', ptemp%id, pblock%id, p
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
end do
|
|
|
|
|
|
|
|
! check neighbors IDs
|
|
|
|
!
|
|
|
|
do i = 1, ndims
|
|
|
|
do j = 1, nsides
|
|
|
|
do k = 1, nfaces
|
|
|
|
ptemp => pblock%neigh(i,j,k)%ptr
|
|
|
|
if (associated(ptemp)) then
|
|
|
|
if (ptemp%id .le. 0 .or. ptemp%id .gt. last_id) then
|
|
|
|
print *, ''
|
|
|
|
print *, ''
|
|
|
|
print *, trim(string)
|
|
|
|
print *, 'wrong neighbor id = ', ptemp%id, pblock%id, i, j, k
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine check_metablock
|
|
|
|
#endif /* DEBUG */
|
2008-11-04 21:00:50 -06:00
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-05 22:16:24 -06:00
|
|
|
!
|
2008-11-04 21:00:50 -06:00
|
|
|
end module
|