2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! This file is part of the AMUN source code, a program to perform
|
|
|
|
!! Newtonian or relativistic magnetohydrodynamical simulations on uniform or
|
|
|
|
!! adaptive mesh.
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2014-01-02 11:52:59 -02:00
|
|
|
!! Copyright (C) 2008-2014 Grzegorz Kowal <grzegorz@amuncode.org>
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! This program is free software: you can redistribute it and/or modify
|
|
|
|
!! it under the terms of the GNU General Public License as published by
|
|
|
|
!! the Free Software Foundation, either version 3 of the License, or
|
|
|
|
!! (at your option) any later version.
|
2008-11-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
|
2012-07-22 12:30:20 -03:00
|
|
|
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-04 21:00:50 -06:00
|
|
|
!!
|
2012-07-22 22:47:25 -03:00
|
|
|
!! module: BLOCKS
|
|
|
|
!!
|
2014-01-15 09:32:58 -02:00
|
|
|
!! This module provides data structures, variables and subroutines to
|
|
|
|
!! construct and dynamically modify the hierarchy of blocks corresponding
|
|
|
|
!! to the simulated mesh geometry.
|
2012-07-22 12:30:20 -03:00
|
|
|
!!
|
|
|
|
!!******************************************************************************
|
2008-11-04 21:00:50 -06:00
|
|
|
!
|
|
|
|
module blocks
|
|
|
|
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! import external subroutines
|
|
|
|
!
|
|
|
|
use timers, only : set_timer, start_timer, stop_timer
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2012-07-22 22:47:25 -03:00
|
|
|
! module variables are not implicit by default
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
2012-07-22 22:47:25 -03:00
|
|
|
implicit none
|
2011-05-05 18:25:24 -03:00
|
|
|
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! timer indices
|
|
|
|
!
|
|
|
|
integer, save :: imi, ima, imu, imp, imq, imr, imd
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! MODULE PARAMETERS:
|
|
|
|
! =================
|
|
|
|
!
|
|
|
|
! ndims - the number of dimensions (2 or 3);
|
|
|
|
! nsides - the number of sides along each direction (2);
|
|
|
|
! nfaces - the number of faces at each side (2 for 2D, 4 for 3D);
|
|
|
|
! nchildren - the number of child blocks for each block (4 for 2D, 8 for 3D);
|
2008-11-04 21:00:50 -06:00
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
integer(kind=4), parameter :: ndims = NDIMS
|
|
|
|
integer(kind=4), parameter :: nsides = 2
|
|
|
|
integer(kind=4), parameter :: nfaces = 2**(ndims - 1)
|
|
|
|
integer(kind=4), parameter :: nchildren = 2**ndims
|
2008-11-04 21:00:50 -06:00
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! MODULE VARIABLES:
|
|
|
|
! ================
|
|
|
|
!
|
|
|
|
! the identification of the last allocated block (always increases)
|
|
|
|
!
|
|
|
|
integer(kind=4), save :: last_id
|
|
|
|
|
|
|
|
! the number of allocated meta and data blocks, and the number of leafs
|
|
|
|
!
|
|
|
|
integer(kind=4), save :: mblocks, dblocks, nleafs
|
|
|
|
|
|
|
|
! the number of variables and fluxes stored in data blocks
|
|
|
|
!
|
|
|
|
integer(kind=4), save :: nvars, nflux
|
|
|
|
|
|
|
|
! the spacial dimensions of allocatable data block arrays
|
|
|
|
!
|
|
|
|
integer(kind=4), save :: nx, ny, nz
|
|
|
|
|
|
|
|
! BLOCK STRUCTURE POINTERS:
|
|
|
|
! ========================
|
|
|
|
!
|
|
|
|
! define pointers to meta, data, and info block structures defined below;
|
|
|
|
! they have to be defined before block structures
|
2009-09-09 16:37:28 -03:00
|
|
|
!
|
|
|
|
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
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! BLOCK STRUCTURES:
|
|
|
|
! ================
|
2009-09-18 17:43:48 -03:00
|
|
|
!
|
2014-01-15 09:32:58 -02:00
|
|
|
! define the META block structure; each process keeps exactly same meta block
|
|
|
|
! structure all the time, so processes can know how the block structure changes
|
|
|
|
! and where to move data blocks;
|
2009-09-09 16:37:28 -03:00
|
|
|
!
|
|
|
|
type block_meta
|
2011-05-05 18:25:24 -03:00
|
|
|
! pointers to the previous and next meta blocks
|
|
|
|
!
|
|
|
|
type(block_meta) , pointer :: prev, next
|
|
|
|
|
|
|
|
! a pointer to the parent meta block
|
|
|
|
!
|
|
|
|
type(block_meta) , pointer :: parent
|
|
|
|
|
2011-05-05 23:36:44 -03:00
|
|
|
! pointers to child meta blocks
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
type(pointer_meta) :: child(nchildren)
|
2011-05-05 18:25:24 -03:00
|
|
|
|
2011-05-05 23:36:44 -03:00
|
|
|
! pointers to neighbor meta blocks
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
type(pointer_meta) :: neigh(ndims,nsides,nfaces)
|
2009-09-09 16:37:28 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
! a pointer to the associated data block
|
|
|
|
!
|
|
|
|
type(block_data) , pointer :: data
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! the identification number (unique for each
|
|
|
|
! block)
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
integer(kind=4) :: id
|
2009-09-09 16:37:28 -03:00
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! the process number to which the meta block
|
|
|
|
! is bounded
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
integer(kind=4) :: cpu
|
|
|
|
|
|
|
|
! the level of refinement
|
|
|
|
!
|
2011-05-05 23:36:44 -03:00
|
|
|
integer(kind=4) :: level
|
2011-05-05 18:25:24 -03:00
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! the number describing the configuration of
|
|
|
|
! the child meta blocks
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
integer(kind=4) :: config
|
|
|
|
|
|
|
|
! the refinement flag, -1, 0, and 1 for
|
2014-01-15 09:32:58 -02:00
|
|
|
! the block marked to be derefined, not
|
|
|
|
! changed, and refined, respectively
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
integer(kind=4) :: refine
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! the position of the block in its siblings
|
|
|
|
! group
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
integer(kind=4) :: pos(ndims)
|
|
|
|
|
|
|
|
! the coordinate of the lower corner of the
|
|
|
|
! block in the effective resolution units
|
|
|
|
!
|
|
|
|
integer(kind=4) :: coord(ndims)
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! the leaf flag, signifying that the block is
|
|
|
|
! the highest block in the local block
|
|
|
|
! structure
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
logical :: leaf
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! the flag indicates that the corresponding
|
|
|
|
! data needs to be updated (e.g. boundaries or
|
|
|
|
! primitive variables), therefore it is
|
|
|
|
! usually .true.
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
2014-01-15 09:32:58 -02:00
|
|
|
logical :: update
|
|
|
|
|
|
|
|
! the block bounds in the coordinate units
|
|
|
|
!
|
|
|
|
real(kind=8) :: xmin, xmax, ymin, ymax, zmin, zmax
|
2011-05-05 18:25:24 -03:00
|
|
|
|
2009-09-09 16:37:28 -03:00
|
|
|
end type block_meta
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! define the DATA block structure; all data blocks are divided between
|
|
|
|
! processes, therefore the same data block cannot be associated with two
|
|
|
|
! different processes, but they can be moved from one process to another;
|
2009-09-09 16:37:28 -03:00
|
|
|
!
|
|
|
|
type block_data
|
2011-05-05 18:25:24 -03:00
|
|
|
! pointers to the previous and next data blocks
|
|
|
|
!
|
|
|
|
type(block_data), pointer :: prev, next
|
|
|
|
|
|
|
|
! a pointer to the associated meta block
|
|
|
|
!
|
|
|
|
type(block_meta), pointer :: meta
|
2009-09-09 16:37:28 -03:00
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! a pointer to the current conserved variable
|
|
|
|
! array
|
2012-07-31 14:14:42 -03:00
|
|
|
!
|
|
|
|
real, dimension(:,:,:,:) , pointer :: u
|
|
|
|
|
|
|
|
! an allocatable arrays to store all conserved
|
2014-01-15 09:32:58 -02:00
|
|
|
! variables (required two for Runge-Kutta
|
|
|
|
! temporal integration methods)
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
2012-07-31 14:14:42 -03:00
|
|
|
real, dimension(:,:,:,:) , allocatable :: u0, u1
|
2011-05-05 18:25:24 -03:00
|
|
|
|
2012-07-27 21:45:35 -03:00
|
|
|
! an allocatable array to store all primitive
|
|
|
|
! variables
|
|
|
|
!
|
|
|
|
real, dimension(:,:,:,:) , allocatable :: q
|
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
! an allocatable array to store all fluxes
|
|
|
|
!
|
|
|
|
real, dimension(:,:,:,:,:), allocatable :: f
|
2009-09-09 16:37:28 -03:00
|
|
|
|
|
|
|
end type block_data
|
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
! define the INFO block structure
|
2009-09-18 17:43:48 -03:00
|
|
|
!
|
|
|
|
type block_info
|
2011-05-05 18:25:24 -03:00
|
|
|
! pointers to the previous and next info blocks
|
|
|
|
!
|
|
|
|
type(block_info) , pointer :: prev, next
|
|
|
|
|
|
|
|
! a pointer to the associated meta block
|
|
|
|
!
|
|
|
|
type(block_meta) , pointer :: block
|
|
|
|
|
|
|
|
! a pointer to the associated neighbor block
|
|
|
|
!
|
|
|
|
type(block_meta) , pointer :: neigh
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! the direction, side and face numbers
|
|
|
|
! indicating the neighbor block orientation
|
|
|
|
! with respect to the block
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
integer(kind=4) :: direction, side, face
|
|
|
|
|
|
|
|
! the level difference between the block and
|
|
|
|
! its neighbor
|
|
|
|
!
|
|
|
|
integer(kind=4) :: level_difference
|
|
|
|
|
2009-09-18 17:43:48 -03:00
|
|
|
end type block_info
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! POINTERS TO THE FIST AND LAST BLOCKS IN THE LISTS:
|
|
|
|
! =================================================
|
|
|
|
!
|
|
|
|
! these pointers construct the lists of meta and data blocks;
|
2009-09-09 16:37:28 -03:00
|
|
|
!
|
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
|
|
|
|
2012-07-22 22:47:25 -03:00
|
|
|
! all variables and subroutines are private by default
|
|
|
|
!
|
|
|
|
private
|
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
! declare public pointers, structures, and variables
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
public :: pointer_meta, pointer_info
|
|
|
|
public :: block_meta, block_data, block_info
|
|
|
|
public :: list_meta, list_data
|
2014-01-15 09:32:58 -02:00
|
|
|
public :: ndims, nsides, nfaces, nchildren
|
|
|
|
|
|
|
|
! declare public subroutines
|
|
|
|
!
|
2012-07-22 22:47:25 -03:00
|
|
|
public :: initialize_blocks, finalize_blocks
|
2014-01-15 11:10:27 -02:00
|
|
|
public :: append_metablock, remove_metablock
|
2014-01-15 11:27:03 -02:00
|
|
|
public :: append_datablock, remove_datablock
|
2014-01-15 11:52:40 -02:00
|
|
|
public :: allocate_metablock, deallocate_metablock
|
2014-01-15 12:04:12 -02:00
|
|
|
public :: allocate_datablock, deallocate_datablock
|
2013-12-23 18:19:13 -02:00
|
|
|
public :: link_blocks, unlink_blocks
|
2014-01-15 12:08:24 -02:00
|
|
|
public :: set_last_id, get_last_id, get_mblocks, get_dblocks, get_nleafs
|
2011-05-06 09:17:01 -03:00
|
|
|
public :: metablock_set_id, metablock_set_cpu, metablock_set_refine &
|
|
|
|
, metablock_set_config, metablock_set_level, metablock_set_position &
|
|
|
|
, metablock_set_coord, metablock_set_bounds, metablock_set_leaf
|
2011-05-05 18:25:24 -03:00
|
|
|
public :: datablock_set_dims
|
|
|
|
public :: refine_block, derefine_block
|
2011-05-06 08:40:21 -03:00
|
|
|
#ifdef DEBUG
|
|
|
|
public :: check_metablock
|
|
|
|
#endif /* DEBUG */
|
2011-05-05 18:25:24 -03:00
|
|
|
|
2014-01-15 09:32:58 -02:00
|
|
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
|
!
|
2008-11-04 21:00:50 -06:00
|
|
|
contains
|
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
!===============================================================================
|
2011-05-05 18:25:24 -03:00
|
|
|
!!
|
2014-01-15 10:33:37 -02:00
|
|
|
!!*** PUBLIC SUBROUTINES *****************************************************
|
2011-05-05 18:25:24 -03:00
|
|
|
!!
|
2014-01-15 10:33:37 -02:00
|
|
|
!===============================================================================
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2012-07-22 22:47:25 -03:00
|
|
|
! subroutine INITIALIZE_BLOCKS:
|
|
|
|
! ----------------------------
|
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
! Subroutine initializes the module structures, pointers and variables.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! verbose - flag determining if the subroutine should be verbose;
|
|
|
|
! iret - return flag of the procedure execution status;
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
subroutine initialize_blocks(verbose, iret)
|
2008-11-07 00:10:09 -06:00
|
|
|
|
2012-07-22 22:47:25 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2008-11-07 00:10:09 -06:00
|
|
|
implicit none
|
2014-01-15 10:33:37 -02:00
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
logical, intent(in) :: verbose
|
|
|
|
integer, intent(inout) :: iret
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! set timer descriptions
|
|
|
|
!
|
|
|
|
call set_timer('blocks:: initialization' , imi)
|
|
|
|
call set_timer('blocks:: meta block allocation' , ima)
|
|
|
|
call set_timer('blocks:: meta block deallocation', imu)
|
|
|
|
call set_timer('blocks:: data block allocation' , imp)
|
|
|
|
call set_timer('blocks:: data block deallocation', imq)
|
|
|
|
call set_timer('blocks:: refine' , imr)
|
|
|
|
call set_timer('blocks:: derefine' , imd)
|
|
|
|
|
|
|
|
! start accounting time for module initialization/finalization
|
|
|
|
!
|
|
|
|
call start_timer(imi)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2014-01-15 10:33:37 -02:00
|
|
|
! reset identification counter
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
last_id = 0
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2011-05-05 23:48:10 -03:00
|
|
|
! reset the number of meta blocks, 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
|
2014-01-15 10:33:37 -02:00
|
|
|
nflux = 1
|
2011-05-05 15:56:38 -03:00
|
|
|
|
|
|
|
! set the initial data block resolution
|
|
|
|
!
|
|
|
|
nx = 1
|
|
|
|
ny = 1
|
|
|
|
nz = 1
|
|
|
|
|
2014-01-15 10:33:37 -02:00
|
|
|
! nullify pointers defining the meta and data lists
|
2008-12-05 14:24:01 -06:00
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
nullify(list_meta)
|
|
|
|
nullify(list_data)
|
|
|
|
nullify(last_meta)
|
|
|
|
nullify(last_data)
|
2011-04-15 00:25:22 +02:00
|
|
|
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! stop accounting time for module initialization/finalization
|
|
|
|
!
|
|
|
|
call stop_timer(imi)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2012-07-22 22:47:25 -03:00
|
|
|
end subroutine initialize_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
|
|
|
!
|
2012-07-22 22:47:25 -03:00
|
|
|
! subroutine FINALIZE_BLOCKS:
|
|
|
|
! --------------------------
|
|
|
|
!
|
|
|
|
! Subroutine iterates over all meta blocks and first deallocates all
|
|
|
|
! associated with them data blocks, and then their metadata structure.
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! iret - return flag of the procedure execution status;
|
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-11-07 00:10:09 -06:00
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
subroutine finalize_blocks(iret)
|
2008-11-07 00:10:09 -06:00
|
|
|
|
2012-07-22 22:47:25 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2008-11-07 00:10:09 -06:00
|
|
|
implicit none
|
|
|
|
|
2014-01-15 10:33:37 -02:00
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer, intent(inout) :: iret
|
|
|
|
|
|
|
|
! local variables
|
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
|
|
|
!
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! start accounting time for module initialization/finalization
|
|
|
|
!
|
|
|
|
call start_timer(imi)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2014-01-15 10:33:37 -02:00
|
|
|
! associate the pointer with the last block on the meta block list
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
pmeta => last_meta
|
2011-05-05 23:48:10 -03:00
|
|
|
|
2014-01-15 10:33:37 -02:00
|
|
|
! iterate until the first block on the list is reached
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
do while(associated(pmeta))
|
2011-05-05 23:48:10 -03:00
|
|
|
|
2014-01-15 10:33:37 -02:00
|
|
|
! deallocate the last meta block
|
2011-05-05 23:48:10 -03:00
|
|
|
!
|
2014-01-15 11:10:27 -02:00
|
|
|
call remove_metablock(pmeta)
|
2011-05-05 23:48:10 -03:00
|
|
|
|
2014-01-15 10:33:37 -02:00
|
|
|
! assign the pointer to the last block on the meta block list
|
2011-05-05 23:48:10 -03:00
|
|
|
!
|
2014-01-15 10:33:37 -02:00
|
|
|
pmeta => last_meta
|
|
|
|
|
|
|
|
end do ! meta blocks
|
2011-05-05 23:48:10 -03:00
|
|
|
|
2014-01-15 10:33:37 -02:00
|
|
|
! nullify pointers defining the meta and data lists
|
|
|
|
!
|
|
|
|
nullify(list_meta)
|
|
|
|
nullify(list_data)
|
|
|
|
nullify(last_meta)
|
|
|
|
nullify(last_data)
|
2010-10-11 22:46:07 -03:00
|
|
|
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! stop accounting time for module initialization/finalization
|
|
|
|
!
|
|
|
|
call stop_timer(imi)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2012-07-22 22:47:25 -03:00
|
|
|
end subroutine finalize_blocks
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2014-01-15 11:10:27 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine APPEND_METABLOCK:
|
|
|
|
! ---------------------------
|
|
|
|
!
|
2014-01-15 11:27:03 -02:00
|
|
|
! Subroutine allocates memory for one meta block, appends it to the meta
|
|
|
|
! block list and returns a pointer associated with it.
|
2014-01-15 11:10:27 -02:00
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
2014-01-15 11:27:03 -02:00
|
|
|
! pmeta - the pointer associated with the newly appended meta block;
|
2014-01-15 11:10:27 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine append_metablock(pmeta)
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(out) :: pmeta
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! allocate memory for the new meta block
|
|
|
|
!
|
|
|
|
call allocate_metablock(pmeta)
|
|
|
|
|
|
|
|
! check if there are any blocks in the meta block list
|
|
|
|
!
|
|
|
|
if (associated(last_meta)) then
|
|
|
|
|
|
|
|
! add the new block to the end of the list
|
|
|
|
!
|
|
|
|
pmeta%prev => last_meta
|
|
|
|
last_meta%next => pmeta
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! there are no blocks in the list, so add this one as the first block
|
|
|
|
!
|
|
|
|
list_meta => pmeta
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
! update the pointer to the last block on the list
|
|
|
|
!
|
|
|
|
last_meta => pmeta
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine append_metablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine REMOVE_METABLOCK:
|
|
|
|
! ---------------------------
|
|
|
|
!
|
2014-01-15 11:27:03 -02:00
|
|
|
! Subroutine removes a meta block associated with the input pointer from
|
|
|
|
! the meta block list, and deallocates space used by it.
|
2014-01-15 11:10:27 -02:00
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
2014-01-15 11:27:03 -02:00
|
|
|
! pmeta - the pointer pointing to the meta block which will be removed;
|
2014-01-15 11:10:27 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine remove_metablock(pmeta)
|
|
|
|
|
|
|
|
! import external procedures
|
|
|
|
!
|
|
|
|
use error , only : print_error
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! check if the pointer is actually associated with any block
|
|
|
|
!
|
|
|
|
if (associated(pmeta)) then
|
|
|
|
|
|
|
|
! if this is the first block in the list, update the list_meta pointer
|
|
|
|
!
|
|
|
|
if (pmeta%id == list_meta%id) list_meta => pmeta%next
|
|
|
|
|
|
|
|
! if this is the last block in the list, update the last_meta pointer
|
|
|
|
!
|
|
|
|
if (pmeta%id == last_meta%id) last_meta => pmeta%prev
|
|
|
|
|
|
|
|
! update the %next and %prev pointers of the previous and next blocks,
|
|
|
|
! respectively
|
|
|
|
!
|
|
|
|
if (associated(pmeta%prev)) pmeta%prev%next => pmeta%next
|
|
|
|
if (associated(pmeta%next)) pmeta%next%prev => pmeta%prev
|
|
|
|
|
|
|
|
! deallocate memory used by the meta block
|
|
|
|
!
|
|
|
|
call deallocate_metablock(pmeta)
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! the argument contains a null pointer, so print an error
|
|
|
|
!
|
|
|
|
call print_error("blocks::remove_metablock" &
|
|
|
|
, "Null pointer argument to meta block!")
|
|
|
|
end if
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine remove_metablock
|
2011-05-05 15:56:38 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2014-01-15 11:27:03 -02:00
|
|
|
! subroutine APPEND_DATABLOCK:
|
|
|
|
! ---------------------------
|
|
|
|
!
|
|
|
|
! Subroutine allocates memory for one data block, appends it to the data
|
|
|
|
! block list and returns a pointer associated with it.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pdata - the pointer associated with the newly appended data block;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine append_datablock(pdata)
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_data), pointer, intent(out) :: pdata
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! allocate memory for the new data block
|
|
|
|
!
|
|
|
|
call allocate_datablock(pdata)
|
|
|
|
|
|
|
|
! check if there are any blocks in the data block list
|
|
|
|
!
|
|
|
|
if (associated(last_data)) then
|
|
|
|
|
|
|
|
! add the new block to the end of the list
|
|
|
|
!
|
|
|
|
pdata%prev => last_data
|
|
|
|
last_data%next => pdata
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! there are no blocks in the list, so add this one as the first block
|
|
|
|
!
|
|
|
|
list_data => pdata
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
! update the pointer to the last block on the list
|
|
|
|
!
|
|
|
|
last_data => pdata
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine append_datablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine REMOVE_DATABLOCK:
|
|
|
|
! ---------------------------
|
|
|
|
!
|
|
|
|
! Subroutine removes a data block associated with the input pointer from
|
|
|
|
! the data block list, and deallocates space used by it.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pdata - the pointer pointing to the data block which will be removed;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine remove_datablock(pdata)
|
|
|
|
|
|
|
|
! import external procedures
|
|
|
|
!
|
|
|
|
use error , only : print_error
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_data), pointer, intent(inout) :: pdata
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! check if the pointer is actually associated with any block
|
|
|
|
!
|
|
|
|
if (associated(pdata)) then
|
|
|
|
|
|
|
|
! check if the data block has associated meta block
|
|
|
|
!
|
|
|
|
if (associated(pdata%meta)) then
|
|
|
|
|
|
|
|
! if this is the first block in the list, update the list_data pointer
|
|
|
|
!
|
|
|
|
if (pdata%meta%id == list_data%meta%id) list_data => pdata%next
|
|
|
|
|
|
|
|
! if this is the last block in the list, update the last_data pointer
|
|
|
|
!
|
|
|
|
if (pdata%meta%id == last_data%meta%id) last_data => pdata%prev
|
|
|
|
|
|
|
|
! update the %next and %prev pointers of the previous and next blocks,
|
|
|
|
! respectively
|
|
|
|
!
|
|
|
|
if (associated(pdata%prev)) pdata%prev%next => pdata%next
|
|
|
|
if (associated(pdata%next)) pdata%next%prev => pdata%prev
|
|
|
|
|
|
|
|
else ! %meta associated
|
|
|
|
|
|
|
|
! there is no meta block associated, so print an error
|
|
|
|
!
|
|
|
|
call print_error("blocks::remove_datablock" &
|
|
|
|
, "No meta block associated with the data block!")
|
|
|
|
|
|
|
|
end if ! %meta associated
|
|
|
|
|
|
|
|
! deallocate the associated data block
|
|
|
|
!
|
|
|
|
call deallocate_datablock(pdata)
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! the argument contains a null pointer, so print an error
|
|
|
|
!
|
|
|
|
call print_error("blocks::remove_datablock" &
|
|
|
|
, "Null pointer argument to data block!")
|
|
|
|
end if
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine remove_datablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
! subroutine ALLOCATE_METABLOCK:
|
|
|
|
! -----------------------------
|
|
|
|
!
|
|
|
|
! Subroutine allocates memory for one meta block, initializes its fields
|
|
|
|
! and returns a pointer associated with it.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pmeta - the pointer associated with the newly allocated meta block;
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine allocate_metablock(pmeta)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
implicit none
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! subroutine arguments
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
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
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
2014-01-15 11:50:48 -02:00
|
|
|
! start accounting time for the meta block allocation
|
2014-01-15 08:42:50 -02:00
|
|
|
!
|
|
|
|
call start_timer(ima)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! allocate the meta block structure for one object
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
allocate(pmeta)
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! nullify fields pointing to previous and next block on the meta block list
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%prev)
|
|
|
|
nullify(pmeta%next)
|
2014-01-15 11:50:48 -02:00
|
|
|
|
|
|
|
! nullify the field pointing to the parent
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%parent)
|
2014-01-15 11:50:48 -02:00
|
|
|
|
|
|
|
! nullify fields pointing to children
|
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do i = 1, nchildren
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%child(i)%ptr)
|
2009-09-10 17:25:28 -03:00
|
|
|
end do
|
2014-01-15 11:50:48 -02:00
|
|
|
|
|
|
|
! nullify fields pointing to neighbors
|
|
|
|
!
|
|
|
|
do i = 1, ndims
|
2009-09-11 21:52:18 -03:00
|
|
|
do j = 1, nsides
|
2014-01-15 11:50:48 -02:00
|
|
|
do k = 1, nfaces
|
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
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! nullify the field pointing to the associated data block
|
|
|
|
!
|
|
|
|
nullify(pmeta%data)
|
|
|
|
|
2009-09-10 17:25:28 -03:00
|
|
|
! set unique ID
|
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
pmeta%id = increase_id()
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! unset the process number, level, the children configuration, refine, leaf,
|
|
|
|
! and update flags
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
pmeta%cpu = -1
|
|
|
|
pmeta%level = -1
|
|
|
|
pmeta%config = -1
|
|
|
|
pmeta%refine = 0
|
|
|
|
pmeta%leaf = .false.
|
|
|
|
pmeta%update = .true.
|
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
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! initialize the effective coordinates
|
2010-10-11 02:10:47 -03:00
|
|
|
!
|
|
|
|
pmeta%coord(:) = 0
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! initialize coordinate bounds of the block
|
2009-09-25 11:41:19 -03:00
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
pmeta%xmin = 0.0d+00
|
|
|
|
pmeta%xmax = 1.0d+00
|
|
|
|
pmeta%ymin = 0.0d+00
|
|
|
|
pmeta%ymax = 1.0d+00
|
|
|
|
pmeta%zmin = 0.0d+00
|
|
|
|
pmeta%zmax = 1.0d+00
|
2009-09-25 11:41:19 -03:00
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! increase the number of allocated meta blocks
|
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
mblocks = mblocks + 1
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
2014-01-15 11:50:48 -02:00
|
|
|
! stop accounting time for the meta block allocation
|
2014-01-15 08:42:50 -02:00
|
|
|
!
|
|
|
|
call stop_timer(ima)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2009-09-10 17:25:28 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine allocate_metablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
! subroutine DEALLOCATE_METABLOCK:
|
|
|
|
! -------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine releases memory used by the meta block associated with
|
|
|
|
! the pointer argument.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pmeta - the pointer associated with the meta block which will be
|
|
|
|
! deallocated;
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
subroutine deallocate_metablock(pmeta)
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! import external procedures
|
|
|
|
!
|
|
|
|
use error , only : print_error
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2009-09-10 16:18:59 -03:00
|
|
|
implicit none
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! subroutine arguments
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
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
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
2014-01-15 11:50:48 -02:00
|
|
|
! start accounting time for the meta block deallocation
|
2014-01-15 08:42:50 -02:00
|
|
|
!
|
|
|
|
call start_timer(imu)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! check if the pointer is actually associated with any block
|
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
if (associated(pmeta)) then
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! decrease the number of leafs
|
|
|
|
!
|
|
|
|
if (pmeta%leaf) nleafs = nleafs - 1
|
|
|
|
|
|
|
|
! decrease the number of allocated meta blocks
|
|
|
|
!
|
|
|
|
mblocks = mblocks - 1
|
|
|
|
|
|
|
|
! nullify fields pointing to previous and next block on the meta block list
|
|
|
|
!
|
|
|
|
nullify(pmeta%prev)
|
|
|
|
nullify(pmeta%next)
|
|
|
|
|
|
|
|
! nullify the field pointing to the parent
|
|
|
|
!
|
|
|
|
nullify(pmeta%parent)
|
|
|
|
|
|
|
|
! nullify fields pointing to children
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do i = 1, nchildren
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%child(i)%ptr)
|
2009-09-10 16:18:59 -03:00
|
|
|
end do
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! nullify fields pointing to neighbors
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
do i = 1, ndims
|
2009-09-11 21:52:18 -03:00
|
|
|
do j = 1, nsides
|
2014-01-15 11:50:48 -02:00
|
|
|
do k = 1, nfaces
|
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
|
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! if there is a data block is associated, remove it
|
2009-09-13 22:58:55 -03:00
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
if (associated(pmeta%data)) call remove_datablock(pmeta%data)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! nullify the field pointing to the associated data block
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
nullify(pmeta%data)
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! release the memory occupied by the block
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2009-09-26 14:27:47 -03:00
|
|
|
deallocate(pmeta)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! nullify the pointer to the deallocated meta block
|
2009-09-13 22:58:55 -03:00
|
|
|
!
|
2014-01-15 11:50:48 -02:00
|
|
|
nullify(pmeta)
|
|
|
|
|
|
|
|
else
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2014-01-15 11:50:48 -02:00
|
|
|
! the argument contains a null pointer, so print an error
|
|
|
|
!
|
|
|
|
call print_error("blocks::deallocate_metablock" &
|
|
|
|
, "Null pointer argument to meta block!")
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
2014-01-15 11:50:48 -02:00
|
|
|
! stop accounting time for the meta block deallocation
|
2014-01-15 08:42:50 -02:00
|
|
|
!
|
|
|
|
call stop_timer(imu)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
2009-09-10 16:18:59 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine deallocate_metablock
|
|
|
|
!
|
2014-01-15 12:04:12 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine ALLOCATE_DATABLOCK:
|
|
|
|
! -----------------------------
|
|
|
|
!
|
|
|
|
! Subroutine allocates memory for one data block, initializes its fields
|
|
|
|
! and returns a pointer associated with it.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pdata - the pointer associated with the newly allocated data block;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine allocate_datablock(pdata)
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_data), pointer, intent(out) :: pdata
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
#ifdef PROFILE
|
|
|
|
! start accounting time for the data block allocation
|
|
|
|
!
|
|
|
|
call start_timer(imp)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
|
|
|
! allocate the block structure
|
|
|
|
!
|
|
|
|
allocate(pdata)
|
|
|
|
|
|
|
|
! nullify field pointing to the previous and next blocks on the data block list
|
|
|
|
!
|
|
|
|
nullify(pdata%prev)
|
|
|
|
nullify(pdata%next)
|
|
|
|
|
|
|
|
! nullify the field pointing to the associate meta block list
|
|
|
|
!
|
|
|
|
nullify(pdata%meta)
|
|
|
|
|
|
|
|
! allocate space for conserved variables
|
|
|
|
!
|
|
|
|
allocate(pdata%u0(nvars,nx,ny,nz), pdata%u1(nvars,nx,ny,nz))
|
|
|
|
|
|
|
|
! allocate space for primitive variables
|
|
|
|
!
|
|
|
|
allocate(pdata%q(nvars,nx,ny,nz))
|
|
|
|
|
|
|
|
! allocate space for numerical fluxes
|
|
|
|
!
|
|
|
|
allocate(pdata%f(ndims,nflux,nx,ny,nz))
|
|
|
|
|
|
|
|
! initiate the conserved variable pointer
|
|
|
|
!
|
|
|
|
pdata%u => pdata%u0
|
|
|
|
|
|
|
|
! increase the number of allocated data blocks
|
|
|
|
!
|
|
|
|
dblocks = dblocks + 1
|
|
|
|
|
|
|
|
#ifdef PROFILE
|
|
|
|
! stop accounting time for the data block allocation
|
|
|
|
!
|
|
|
|
call stop_timer(imp)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine allocate_datablock
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine DEALLOCATE_DATABLOCK:
|
|
|
|
! -------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine releases memory used by the data block associated with
|
|
|
|
! the pointer argument.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pdata - the pointer associated with the data block which will be
|
|
|
|
! deallocated;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine deallocate_datablock(pdata)
|
|
|
|
|
|
|
|
! import external procedures
|
|
|
|
!
|
|
|
|
use error , only : print_error
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_data), pointer, intent(inout) :: pdata
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
#ifdef PROFILE
|
|
|
|
! start accounting time for the data block deallocation
|
|
|
|
!
|
|
|
|
call start_timer(imq)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
|
|
|
! check if the pointer is actually associated with any block
|
|
|
|
!
|
|
|
|
if (associated(pdata)) then
|
|
|
|
|
|
|
|
! decrease the number of allocated data blocks
|
|
|
|
!
|
|
|
|
dblocks = dblocks - 1
|
|
|
|
|
|
|
|
! nullify field pointing to the previous and next blocks on the data block list
|
|
|
|
!
|
|
|
|
nullify(pdata%prev)
|
|
|
|
nullify(pdata%next)
|
|
|
|
|
|
|
|
! nullify the field pointing to the associate meta block list
|
|
|
|
!
|
|
|
|
nullify(pdata%meta)
|
|
|
|
|
|
|
|
! nullify pointer to the current conserved variable array
|
|
|
|
!
|
|
|
|
nullify(pdata%u)
|
|
|
|
|
|
|
|
! deallocate conserved variables
|
|
|
|
!
|
|
|
|
if (allocated(pdata%u0)) deallocate(pdata%u0)
|
|
|
|
if (allocated(pdata%u1)) deallocate(pdata%u1)
|
|
|
|
|
|
|
|
! deallocate primitive variables
|
|
|
|
!
|
|
|
|
if (allocated(pdata%q )) deallocate(pdata%q )
|
|
|
|
|
|
|
|
! deallocate numerical fluxes
|
|
|
|
!
|
|
|
|
if (allocated(pdata%f )) deallocate(pdata%f )
|
|
|
|
|
|
|
|
! release the memory occupied by the block
|
|
|
|
!
|
|
|
|
deallocate(pdata)
|
|
|
|
|
|
|
|
! nullify the pointer to the deallocated meta block
|
|
|
|
!
|
|
|
|
nullify(pdata)
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! the argument contains a null pointer, so print an error
|
|
|
|
!
|
|
|
|
call print_error("blocks::deallocate_datablock" &
|
|
|
|
, "Null pointer argument to data block!")
|
|
|
|
|
|
|
|
end if ! pdata associated with a data block
|
|
|
|
|
|
|
|
#ifdef PROFILE
|
|
|
|
! stop accounting time for the data block deallocation
|
|
|
|
!
|
|
|
|
call stop_timer(imq)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine deallocate_datablock
|
|
|
|
!
|
2014-01-15 12:08:24 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine LINK_BLOCKS:
|
|
|
|
! ----------------------
|
|
|
|
!
|
|
|
|
! Subroutine links meta and data blocks.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pmeta - the meta block pointer;
|
|
|
|
! pdata - the data block pointer;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine link_blocks(pmeta, pdata)
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
type(block_data), pointer, intent(inout) :: pdata
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! associate the corresponging pointers
|
|
|
|
!
|
|
|
|
pmeta%data => pdata
|
|
|
|
pdata%meta => pmeta
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine link_blocks
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine UNLINK_BLOCKS:
|
|
|
|
! ------------------------
|
|
|
|
!
|
|
|
|
! Subroutine unlinks meta and data blocks.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pmeta - the meta block pointer;
|
|
|
|
! pdata - the data block pointer;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine unlink_blocks(pmeta, pdata)
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
type(block_data), pointer, intent(inout) :: pdata
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! nullify the corresponging pointers
|
|
|
|
!
|
|
|
|
nullify(pmeta%data)
|
|
|
|
nullify(pdata%meta)
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine unlink_blocks
|
|
|
|
!
|
2014-01-15 14:13:31 -02:00
|
|
|
!===============================================================================
|
|
|
|
!!
|
|
|
|
!!*** PRIVATE SUBROUTINES ****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine INSERT_METABLOCK_AFTER:
|
|
|
|
! ---------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine allocates memory for one meta block, inserts it to the meta
|
|
|
|
! block list after the provided pointer and returns a pointer associated
|
|
|
|
! with it.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pmeta - the pointer associated with the newly appended meta block;
|
|
|
|
! pprev - the pointer after which the new block has to be inserted;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine insert_metablock_after(pprev, pmeta)
|
|
|
|
|
|
|
|
! import external procedures
|
|
|
|
!
|
|
|
|
use error , only : print_error
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(in) :: pprev
|
|
|
|
type(block_meta), pointer, intent(out) :: pmeta
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! allocate memory for the new meta block
|
|
|
|
!
|
|
|
|
call allocate_metablock(pmeta)
|
|
|
|
|
|
|
|
! if pprev is associated, insert the new block after it
|
|
|
|
!
|
|
|
|
if (associated(pprev)) then
|
|
|
|
|
|
|
|
! associate the %prev and %next pointers
|
|
|
|
!
|
|
|
|
pmeta%prev => pprev
|
|
|
|
pmeta%next => pprev%next
|
|
|
|
|
|
|
|
! update the pointer of the next and previous blocks
|
|
|
|
!
|
|
|
|
if (associated(pprev%next)) pprev%next%prev => pmeta
|
|
|
|
pprev%next => pmeta
|
|
|
|
|
|
|
|
! check if last_meta is associated
|
|
|
|
!
|
|
|
|
if (associated(last_meta)) then
|
|
|
|
|
|
|
|
! update the last_meta pointer if necessary
|
|
|
|
!
|
|
|
|
if (pprev%id == last_meta%id) last_meta => pmeta
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! strange situation, pprev is associated, but last_meta not
|
|
|
|
!
|
|
|
|
call print_error("blocks::intert_metablock_after" &
|
|
|
|
, "Argument pprev is associated but last_meta is not!")
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! if pprev is null and list_meta is associated, there is something wrong
|
|
|
|
!
|
|
|
|
if (associated(list_meta)) then
|
|
|
|
|
|
|
|
! strange situation, pprev is null but list_meta is associated
|
|
|
|
!
|
|
|
|
call print_error("blocks::intert_metablock_after" &
|
|
|
|
, "Argument pprev is null but list_meta is associated!")
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! pprev and list_meta are nulls, so add the first block to the list by
|
|
|
|
! associating list_meta and last_meta
|
|
|
|
!
|
|
|
|
list_meta => pmeta
|
|
|
|
last_meta => pmeta
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine insert_metablock_after
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine INSERT_METABLOCK_BEFORE:
|
|
|
|
! ----------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine allocates memory for one meta block, inserts it to the meta
|
|
|
|
! block list before the provided pointer and returns a pointer associated
|
|
|
|
! with it.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! pmeta - the pointer associated with the newly appended meta block;
|
|
|
|
! pnext - the pointer before which the new block has to be inserted;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine insert_metablock_before(pnext, pmeta)
|
|
|
|
|
|
|
|
! import external procedures
|
|
|
|
!
|
|
|
|
use error , only : print_error
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(in) :: pnext
|
|
|
|
type(block_meta), pointer, intent(out) :: pmeta
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! allocate memory for the new meta block
|
|
|
|
!
|
|
|
|
call allocate_metablock(pmeta)
|
|
|
|
|
|
|
|
! if pnext is associated, insert the new block before it
|
|
|
|
!
|
|
|
|
if (associated(pnext)) then
|
|
|
|
|
|
|
|
! associate the %prev and %next pointers
|
|
|
|
!
|
|
|
|
pmeta%prev => pnext%prev
|
|
|
|
pmeta%next => pnext
|
|
|
|
|
|
|
|
! update the pointer of the next and previous blocks
|
|
|
|
!
|
|
|
|
if (associated(pnext%prev)) pnext%prev%next => pmeta
|
|
|
|
pnext%prev => pmeta
|
|
|
|
|
|
|
|
! check if list_meta is associated
|
|
|
|
!
|
|
|
|
if (associated(list_meta)) then
|
|
|
|
|
|
|
|
! update the list_meta pointer if necessary
|
|
|
|
!
|
|
|
|
if (pnext%id == list_meta%id) list_meta => pmeta
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! strange situation, pnext is associated, but list_meta not
|
|
|
|
!
|
|
|
|
call print_error("blocks::intert_metablock_before" &
|
|
|
|
, "Argument pnext is associated but list_meta is not!")
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! if pnext is null and last_meta is associated, there is something wrong
|
|
|
|
!
|
|
|
|
if (associated(last_meta)) then
|
|
|
|
|
|
|
|
! strange situation, pnext is null but last_meta is associated
|
|
|
|
!
|
|
|
|
call print_error("blocks::intert_metablock_before" &
|
|
|
|
, "Argument pnext is null but last_meta is associated!")
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! pnext and last_meta are nulls, so add the first block to the list by
|
|
|
|
! associating list_meta and last_meta
|
|
|
|
!
|
|
|
|
list_meta => pmeta
|
|
|
|
last_meta => pmeta
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine insert_metablock_before
|
|
|
|
!
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
!!==============================================================================
|
|
|
|
!!
|
|
|
|
!! TOOL SUBROUTINES
|
|
|
|
!!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 23:36:44 -03:00
|
|
|
! increase_id: function increases the last identification by 1 and returns its
|
2011-05-05 18:25:24 -03:00
|
|
|
! value
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
function increase_id()
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
implicit none
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
! return variable
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
integer(kind=4) :: increase_id
|
2011-04-30 12:11:57 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
! increase ID by 1
|
|
|
|
!
|
|
|
|
last_id = last_id + 1
|
2009-09-10 16:18:59 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
! return ID
|
2009-09-10 16:18:59 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
increase_id = last_id
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
return
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function increase_id
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 23:36:44 -03:00
|
|
|
! set_last_id: subroutine sets the last identification value
|
2011-05-05 18:25:24 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine set_last_id(id)
|
2009-09-13 22:58:55 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
use error, only : print_error
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input argument
|
|
|
|
!
|
|
|
|
integer(kind=4), intent(in) :: id
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
if (last_id .gt. id) then
|
|
|
|
call print_error("blocks::set_last_id" &
|
|
|
|
, "New last_id must be larger than old one!")
|
|
|
|
else
|
|
|
|
last_id = id
|
|
|
|
end if
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine set_last_id
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2008-12-18 23:47:58 -06:00
|
|
|
!===============================================================================
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2011-05-05 23:36:44 -03:00
|
|
|
! get_last_id: function returns the last identification value
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
function get_last_id()
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
! return variable
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
integer(kind=4) :: get_last_id
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
get_last_id = last_id
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
return
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function get_last_id
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! get_mblocks: function returns the number of meta blocks
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
function get_mblocks()
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! return variable
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
integer(kind=4) :: get_mblocks
|
2009-09-21 19:02:05 -03:00
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
get_mblocks = mblocks
|
|
|
|
|
|
|
|
return
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function get_mblocks
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
! get_dblocks: function returns the number of data blocks
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
function get_dblocks()
|
2009-09-10 17:25:28 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
! return variable
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
integer(kind=4) :: get_dblocks
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
get_dblocks = dblocks
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
return
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function get_dblocks
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! get_nleafs: function returns the number of leafs
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
function get_nleafs()
|
2009-09-10 17:25:28 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! return variable
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
integer(kind=4) :: get_nleafs
|
2009-09-21 19:02:05 -03:00
|
|
|
!
|
2009-09-10 17:25:28 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
get_nleafs = nleafs
|
|
|
|
|
|
|
|
return
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end function get_nleafs
|
2009-09-10 17:25:28 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-06 09:17:01 -03:00
|
|
|
! metablock_set_id: subroutine sets the identification value
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine metablock_set_id(pmeta, id)
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
integer(kind=4) , intent(in) :: id
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set the id field
|
|
|
|
!
|
|
|
|
pmeta%id = id
|
|
|
|
|
|
|
|
! check if the id is larger then last_id, if so reset last_id to id
|
|
|
|
!
|
|
|
|
if (last_id .lt. id) last_id = id
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_set_id
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! metablock_set_cpu: subroutine sets the cpu number
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine metablock_set_cpu(pmeta, cpu)
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
integer(kind=4) , intent(in) :: cpu
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! set the cpu field
|
|
|
|
!
|
|
|
|
pmeta%cpu = cpu
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_set_cpu
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! metablock_set_refine: subroutine sets the refine flag
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine metablock_set_refine(pmeta, refine)
|
|
|
|
|
|
|
|
use error, only : print_error
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
type(block_meta), pointer, intent(inout) :: pmeta
|
|
|
|
integer(kind=4) , intent(in) :: refine
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! check if the refine value is correct
|
|
|
|
!
|
|
|
|
if (abs(refine) .gt. 1) then
|
|
|
|
|
|
|
|
! print error about wrong refine flag
|
|
|
|
!
|
|
|
|
call print_error("blocks::metablock_set_refine" &
|
|
|
|
, "New refine flag is incorrect!")
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! set the refine field
|
|
|
|
!
|
|
|
|
pmeta%refine = refine
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_set_refine
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! metablock_set_leaf: subroutine marks the block as a leaf
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
subroutine metablock_set_leaf(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
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2009-09-10 18:15:30 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
end subroutine metablock_set_leaf
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-06 09:17:01 -03:00
|
|
|
! metablock_unset_leaf: subroutine unmarks the block as a leaf
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
subroutine metablock_unset_leaf(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
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2009-09-10 18:15:30 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
end subroutine metablock_unset_leaf
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-06 09:17:01 -03:00
|
|
|
! metablock_set_config: subroutine sets the configuration flag
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
subroutine metablock_set_config(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
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2009-09-10 18:15:30 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
end subroutine metablock_set_config
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
! metablock_set_level: subroutine sets the level of data block
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
subroutine metablock_set_level(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
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2009-09-10 18:15:30 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
end subroutine metablock_set_level
|
2009-09-10 18:15:30 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
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
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 23:36:44 -03:00
|
|
|
! set the coordinates
|
2010-10-11 02:10:47 -03:00
|
|
|
!
|
|
|
|
pmeta%coord(1) = px
|
|
|
|
pmeta%coord(2) = py
|
|
|
|
#if NDIMS == 3
|
|
|
|
pmeta%coord(3) = pz
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine metablock_set_coord
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
! metablock_set_bounds: subroutine sets the bounds of data block
|
2009-09-10 17:46:36 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
subroutine metablock_set_bounds(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
|
2011-05-05 23:36:44 -03:00
|
|
|
real , intent(in) :: xmin, xmax
|
|
|
|
real , intent(in) :: ymin, ymax
|
|
|
|
real , intent(in) :: 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
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2009-09-10 17:46:36 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
end subroutine metablock_set_bounds
|
|
|
|
!
|
2013-12-23 17:40:49 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
! datablock_set_dims: subroutine sets the number of variables and dimensions
|
|
|
|
! for arrays allocated in data blocks
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine datablock_set_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
|
2012-07-27 13:14:48 -03:00
|
|
|
#if NDIMS == 3
|
2011-05-05 18:25:24 -03:00
|
|
|
nz = nk
|
2012-07-27 13:14:48 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine datablock_set_dims
|
|
|
|
!
|
|
|
|
!!==============================================================================
|
|
|
|
!!
|
|
|
|
!! REFINEMENT/DEREFINEMENT SUBROUTINES
|
|
|
|
!!
|
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
|
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
integer, dimension(nchildren) :: config, order
|
2009-09-13 22:58:55 -03:00
|
|
|
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
|
|
|
!
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! start accounting time for the block refinement
|
|
|
|
!
|
|
|
|
call start_timer(imr)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
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
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
call metablock_unset_leaf(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
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do p = 1, nchildren
|
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
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
call metablock_set_leaf(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
|
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do p = 1, nchildren
|
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
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do p = 2, nchildren
|
2009-09-11 21:52:18 -03:00
|
|
|
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
|
2013-12-27 17:58:40 -02:00
|
|
|
pfirst => pblock%child(order( 1))%ptr
|
|
|
|
plast => pblock%child(order(nchildren))%ptr
|
2009-09-11 21:52:18 -03:00
|
|
|
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
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do p = 1, nchildren
|
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
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
call metablock_set_bounds(pchild, xmn, xmx, ymn, ymx, zmn, zmx)
|
2009-09-25 11:41:19 -03:00
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! allocate data blocks if necessary
|
|
|
|
!
|
|
|
|
if (falloc_data) then
|
|
|
|
|
|
|
|
! iterate over all children and allocate data blocks
|
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do p = 1, nchildren
|
2009-09-25 11:41:19 -03:00
|
|
|
|
|
|
|
! 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
|
|
|
|
!
|
2013-12-23 18:19:13 -02:00
|
|
|
call link_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
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do p = 2, nchildren
|
2011-05-05 23:36:44 -03:00
|
|
|
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
|
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
|
|
|
! 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
|
|
|
|
2013-12-27 17:58:40 -02:00
|
|
|
pfirst => pblock%child(order( 1))%ptr
|
|
|
|
plast => pblock%child(order(nchildren))%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 23:36:44 -03:00
|
|
|
call print_error("blocks::refine_block" &
|
|
|
|
, "Input pointer is not associated! Terminating!")
|
|
|
|
|
2009-09-21 19:02:05 -03:00
|
|
|
end if
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! stop accounting time for the block refinement
|
|
|
|
!
|
|
|
|
call stop_timer(imr)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
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
|
|
|
!
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! start accounting time for the block derefinement
|
|
|
|
!
|
|
|
|
call start_timer(imd)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
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
|
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do p = 1, nchildren
|
2011-05-05 18:25:24 -03:00
|
|
|
call metablock_unset_leaf(pblock%child(p)%ptr)
|
2014-01-15 11:10:27 -02:00
|
|
|
call remove_metablock(pblock%child(p)%ptr)
|
2009-09-29 17:19:26 -03:00
|
|
|
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
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
call metablock_set_leaf(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
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2014-01-15 08:42:50 -02:00
|
|
|
#ifdef PROFILE
|
|
|
|
! stop accounting time for the block derefinement
|
|
|
|
!
|
|
|
|
call stop_timer(imd)
|
|
|
|
#endif /* PROFILE */
|
|
|
|
|
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
|
2011-05-05 18:25:24 -03:00
|
|
|
!!==============================================================================
|
|
|
|
!!
|
|
|
|
!! DEBUG SUBROUTINES
|
|
|
|
!!
|
2011-04-21 14:13:25 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! 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
|
|
|
|
!
|
2013-12-27 17:58:40 -02:00
|
|
|
do p = 1, nchildren
|
2011-04-21 14:13:25 -03:00
|
|
|
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
|
2011-05-05 23:41:42 -03:00
|
|
|
|
2011-04-21 14:13:25 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
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
|