!!****************************************************************************** !! !! This file is part of the AMUN source code, a program to perform !! Newtonian or relativistic magnetohydrodynamical simulations on uniform or !! adaptive mesh. !! !! Copyright (C) 2008-2023 Grzegorz Kowal !! !! This program is free software: you can redistribute it and/or modify !! it under the terms of the GNU General Public License as published by !! the Free Software Foundation, either version 3 of the License, or !! (at your option) any later version. !! !! This program is distributed in the hope that it will be useful, !! but WITHOUT ANY WARRANTY; without even the implied warranty of !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the !! GNU General Public License for more details. !! !! You should have received a copy of the GNU General Public License !! along with this program. If not, see . !! !!****************************************************************************** !! !! module: BLOCKS !! !! This module provides data structures, variables and subroutines to !! construct and dynamically modify the hierarchy of blocks corresponding !! to the simulated mesh geometry. !! !!****************************************************************************** ! module blocks implicit none ! MODULE PARAMETERS: ! ================= ! ! ndims - the number of dimensions (2 or 3); ! nsides - the number of sides along each direction (2); ! nchildren - the number of child blocks for each block (4 for 2D, 8 for 3D); ! integer(kind=4), parameter :: ndims = NDIMS integer(kind=4), parameter :: nsides = 2 integer(kind=4), parameter :: nchildren = 2**ndims ! MODULE VARIABLES: ! ================ ! ! the identification of the last allocated block (always increases) ! integer(kind=4), save :: last_id = 0 ! the number of allocated meta and data blocks (inserted in the lists), ! and the number of leafs ! integer(kind=4), save :: mblocks = 0, dblocks = 0, nleafs = 0 ! the number of variables, fluxes, and registers stored in data blocks ! integer(kind=4), save :: nvars = 1, nflux = 1, nregs = 1 ! the spacial dimensions of allocatable data block arrays ! integer(kind=4), save :: nx = 1, ny = 1, nz = 1 ! BLOCK STRUCTURE POINTERS: ! ======================== ! ! define pointers to meta, data, and info block structures defined below; ! they have to be defined before block structures ! type pointer_meta type(block_meta), pointer :: ptr end type pointer_meta type pointer_data type(block_data), pointer :: ptr end type pointer_data type pointer_info type(block_info), pointer :: ptr end type pointer_info ! BLOCK STRUCTURES: ! ================ ! ! 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; ! type block_meta ! 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 ! pointers to child meta blocks ! type(pointer_meta) :: child(nchildren) #if NDIMS == 2 ! pointers to edge neighbor meta blocks with ! indices: ! 1 - the X corner coordinate ! 2 - the Y corner coordinate ! 3 - the direction of the edge from the corner ! with above coordinates ! and dimensions [1:2,1:2,1:2] ! type(pointer_meta) :: edges(nsides,nsides,ndims) ! pointers to corner neighbor meta blocks with ! indices: ! 1 - the X corner coordinate ! 2 - the Y corner coordinate ! and dimensions [1:2,1:2] ! type(pointer_meta) :: corners(nsides,nsides) #endif /* NDIMS == 2 */ #if NDIMS == 3 ! pointers to face neighbor meta blocks with ! indices: ! 1 - the X corner coordinate ! 2 - the Y corner coordinate ! 3 - the Z corner coordinate ! 4 - the direction of the face normal vector ! from the corner with above coordinates ! and dimensions [1:2,1:2,1:2,1:3] ! type(pointer_meta) :: faces(nsides,nsides,nsides,ndims) ! pointers to edge neighbor meta blocks with ! indices: ! 1 - the X corner coordinate ! 2 - the Y corner coordinate ! 3 - the Z corner coordinate ! 4 - the direction of the edge from the corner ! with above coordinates ! and dimensions [1:2,1:2,1:2,1:3] ! type(pointer_meta) :: edges(nsides,nsides,nsides,ndims) ! pointers to corner neighbor meta blocks with ! indices: ! 1 - the X corner coordinate ! 2 - the Y corner coordinate ! 3 - the Z corner coordinate ! and dimensions [1:2,1:2,1:2] ! type(pointer_meta) :: corners(nsides,nsides,nsides) #endif /* NDIMS == 3 */ ! a pointer to the associated data block ! type(block_data) , pointer :: data ! the identification number (unique for each ! block) ! integer(kind=4) :: id ! the process number to which the meta block ! is bounded ! integer(kind=4) :: process ! the level of refinement ! integer(kind=4) :: level ! the number describing the configuration of ! the child meta blocks ! integer(kind=4) :: conf ! the refinement flag, -1, 0, and 1 for ! the block marked to be derefined, not ! changed, and refined, respectively ! integer(kind=4) :: refine ! the block position in its parent ! integer(kind=4) :: pos(ndims) ! the block global coordinates at its level ! integer(kind=4) :: coords(ndims) ! the leaf flag, signifying that the block is ! the highest block in the local block ! structure ! logical :: leaf ! the flag indicates that the corresponding ! block needs to be updated since it was ! refined or derefined ! logical :: update ! the flag indicates that the boundary of ! the block have been updated, so it needs ! the update of the primitive variables ! logical :: boundary ! the block bounds in the coordinate units ! real(kind=8) :: xmin, xmax, ymin, ymax, zmin, zmax end type block_meta ! define the LEAF structure; this is a simple structure used to keep the leaf ! list; it contains the pointer to the next leaf block and the pointer to ! the meta block which is leaf; ! type block_leaf ! pointer to the next leaf block ! type(block_leaf), pointer :: next ! pointer to the leaf meta block ! type(block_meta), pointer :: meta end type block_leaf ! 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; ! type block_data ! 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 ! a pointer to the current conserved variable ! array ! real(kind=8), dimension(:,:,:,:) , pointer :: u ! an allocatable arrays to store all conserved ! variables (requires two additional arrays ! for embedded low-storage explicit optimal ! high-order SSPRK integration methods) ! real(kind=8), dimension(:,:,:,:,:), allocatable :: uu ! an allocatable array to store all primitive ! variables and ! real(kind=8), dimension(:,:,:,:) , allocatable :: q, du ! allocatable arrays to store interface fluxes ! real(kind=8), dimension(:,:,:,:) , allocatable :: fx, fy, fz ! flag indicating if all block cells ! are physical ! logical :: physical end type block_data ! define the INFO block structure ! type block_info ! 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 :: meta ! a pointer to the associated neighbor block ! type(block_meta) , pointer :: neigh ! the direction along which the neighbor ! is located ! integer(kind=4) :: direction ! the boundary region location ! integer(kind=4) :: location(NDIMS) end type block_info ! POINTERS TO THE FIST AND LAST BLOCKS IN THE LISTS: ! ================================================= ! ! these pointers construct the lists of meta and data blocks, and the list of ! leaf blocks; ! type(block_meta), pointer, save :: list_meta, last_meta type(block_data), pointer, save :: list_data, last_data type(block_leaf), pointer, save :: list_leaf ! THE VECTOR OF DATA BLOCK POINTERS: ! ================================= ! type(pointer_data), dimension(:), allocatable :: data_blocks ! all variables and subroutines are private by default ! private ! declare public pointers, structures, and variables ! public :: pointer_meta, pointer_data, pointer_info public :: block_meta, block_data, block_info, block_leaf public :: list_meta, list_data, list_leaf public :: ndims, nsides, nchildren, nregs ! declare public subroutines ! public :: initialize_blocks, finalize_blocks public :: append_metablock, remove_metablock public :: append_datablock, remove_datablock public :: allocate_metablock, deallocate_metablock public :: allocate_datablock, deallocate_datablock public :: link_blocks, unlink_blocks public :: refine_block, derefine_block public :: set_last_id, get_last_id, get_mblocks, get_dblocks, get_nleafs public :: set_blocks_update public :: change_blocks_process public :: set_neighbors_refine, set_neighbors_update public :: metablock_set_id, metablock_set_process, metablock_set_level public :: metablock_set_configuration, metablock_set_refinement public :: metablock_set_position, metablock_set_coordinates public :: metablock_set_bounds, metablock_set_leaf, metablock_unset_leaf public :: build_leaf_list, build_datablock_list public :: data_blocks #ifdef DEBUG public :: check_neighbors #endif /* DEBUG */ !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! contains ! !=============================================================================== !! !!*** PUBLIC SUBROUTINES ***************************************************** !! !=============================================================================== ! !=============================================================================== ! ! subroutine INITIALIZE_BLOCKS: ! ---------------------------- ! ! Subroutine initializes the module structures, pointers and variables. ! ! Arguments: ! ! nv - the number of variables stored in the block; ! nf - the number of fluxes stored in the block; ! nc - the number of cells per block dimension; ! nr - the number of registers; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine initialize_blocks(nv, nf, nc, nr, status) implicit none integer, intent(in) :: nv, nf, nc, nr integer, intent(out) :: status !------------------------------------------------------------------------------- ! status = 0 ! reset identification counter ! last_id = 0 ! reset the number of meta blocks, data blocks, and leafs ! mblocks = 0 dblocks = 0 nleafs = 0 ! set the initial number of variables, fluxes, and registers ! nvars = nv nflux = nf nregs = nr ! set the initial data block resolution ! nx = nc ny = nc #if NDIMS == 3 nz = nc #endif /* NDIMS == 3 */ ! nullify pointers defining the meta and data lists ! nullify(list_meta) nullify(list_data) nullify(last_meta) nullify(last_data) nullify(list_leaf) !------------------------------------------------------------------------------- ! end subroutine initialize_blocks ! !=============================================================================== ! ! subroutine FINALIZE_BLOCKS: ! -------------------------- ! ! Subroutine iterates over all meta blocks and first deallocates all ! associated with them data blocks, and then their metadata structure. ! ! Arguments: ! ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine finalize_blocks(status) use iso_fortran_env, only : error_unit implicit none integer, intent(out) :: status type(block_meta), pointer :: pmeta character(len=*), parameter :: loc = 'BLOCKS::finalize_blocks()' !------------------------------------------------------------------------------- ! status = 0 ! deallocate the vector of data blocks ! if (allocated(data_blocks)) then deallocate(data_blocks, stat=status) if (status /= 0) & write(error_unit,"('[',a,']: ',a)") trim(loc), & "Could not deallocate the vector of data blocks!" end if ! remove all metablocks ! pmeta => last_meta do while(associated(pmeta)) call remove_metablock(pmeta, status) pmeta => last_meta end do ! meta blocks ! wipe the leaf list out ! call wipe_leaf_list(status) ! nullify pointers defining the meta and data lists ! nullify(list_meta) nullify(list_data) nullify(last_meta) nullify(last_data) nullify(list_leaf) !------------------------------------------------------------------------------- ! end subroutine finalize_blocks ! !=============================================================================== ! ! subroutine APPEND_METABLOCK: ! --------------------------- ! ! Subroutine allocates memory for one meta block, appends it to the meta ! block list and returns a pointer associated with it. ! ! Arguments: ! ! pmeta - the pointer associated with the newly appended meta block; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine append_metablock(pmeta, status) implicit none type(block_meta), pointer, intent(out) :: pmeta integer , intent(out) :: status ! !------------------------------------------------------------------------------- ! ! allocate memory for the new meta block ! call allocate_metablock(pmeta, status) ! quit if block couldn't be allocated ! if (status /= 0) return ! 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 ! increase the number of allocated meta blocks stored in the meta block list ! mblocks = mblocks + 1 !------------------------------------------------------------------------------- ! end subroutine append_metablock ! !=============================================================================== ! ! subroutine REMOVE_METABLOCK: ! --------------------------- ! ! Subroutine removes a meta block associated with the input pointer from ! the meta block list, and deallocates space used by it. ! ! Arguments: ! ! pmeta - the pointer pointing to the meta block which will be removed; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine remove_metablock(pmeta, status) use iso_fortran_env, only : error_unit implicit none type(block_meta), pointer, intent(inout) :: pmeta integer , intent(out) :: status character(len=*), parameter :: loc = 'BLOCKS::remove_metablock()' !------------------------------------------------------------------------------- ! ! 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 ! set this block to be not a leaf ! call metablock_unset_leaf(pmeta) ! decrease the number of allocated meta blocks stored in the meta block list ! mblocks = mblocks - 1 ! deallocate memory used by the meta block ! call deallocate_metablock(pmeta, status) else ! the argument contains a null pointer, so print an error ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Null pointer argument to meta block!" status = 1 end if !------------------------------------------------------------------------------- ! end subroutine remove_metablock ! !=============================================================================== ! ! 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; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine append_datablock(pdata, status) implicit none type(block_data), pointer, intent(out) :: pdata integer , intent(out) :: status !------------------------------------------------------------------------------- ! ! allocate memory for the new data block ! call allocate_datablock(pdata, status) ! quit if block couldn't be allocated ! if (status /= 0) return ! 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 ! increase the number of data blocks in the list ! dblocks = dblocks + 1 !------------------------------------------------------------------------------- ! 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; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine remove_datablock(pdata, status) use iso_fortran_env, only : error_unit implicit none type(block_data), pointer, intent(inout) :: pdata integer , intent(out) :: status character(len=*), parameter :: loc = 'BLOCKS::remove_datablock()' !------------------------------------------------------------------------------- ! ! 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 ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "No meta block associated with the data block!" status = 1 return end if ! %meta associated ! decrease the number of allocated data blocks in the list ! dblocks = dblocks - 1 ! deallocate the associated data block ! call deallocate_datablock(pdata, status) else ! the argument contains a null pointer, so print an error ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Null pointer argument to data block!" status = 1 end if !------------------------------------------------------------------------------- ! end subroutine remove_datablock ! !=============================================================================== ! ! 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; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine allocate_metablock(pmeta, status) use iso_fortran_env, only : error_unit implicit none type(block_meta), pointer, intent(out) :: pmeta integer , intent(out) :: status integer :: n, i, j #if NDIMS == 3 integer :: k #endif /* NDIMS == 3 */ character(len=*), parameter :: loc = 'BLOCKS::allocate_metablock()' !------------------------------------------------------------------------------- ! ! allocate the meta block structure for one object ! allocate(pmeta, stat = status) ! check if allocation succeeded ! if (status == 0) then ! 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 ! do i = 1, nchildren nullify(pmeta%child(i)%ptr) end do ! nullify fields pointing to face, edge, and corner neighbors ! #if NDIMS == 2 do i = 1, nsides do j = 1, nsides do n = 1, ndims nullify(pmeta%edges(i,j,n)%ptr) end do ! ndims nullify(pmeta%corners(i,j)%ptr) end do ! nsides end do ! nsides #endif /* NDIMS == 2 */ #if NDIMS == 3 do i = 1, nsides do j = 1, nsides do k = 1, nsides do n = 1, ndims nullify(pmeta%faces(i,j,k,n)%ptr) nullify(pmeta%edges(i,j,k,n)%ptr) end do ! ndims nullify(pmeta%corners(i,j,k)%ptr) end do ! nsides end do ! nsides end do ! nsides #endif /* NDIMS == 3 */ ! nullify the field pointing to the associated data block ! nullify(pmeta%data) ! set unique ID ! pmeta%id = increase_id() ! unset the process number, level, the children configuration, refine, leaf, ! and update flags ! pmeta%process = -1 pmeta%level = -1 pmeta%conf = -1 pmeta%refine = 0 pmeta%leaf = .false. pmeta%update = .true. pmeta%boundary = .false. ! initialize the position in the parent block ! pmeta%pos(:) = -1 ! initialize the block coordinates in the current level ! pmeta%coords(:) = 0 ! initialize coordinate bounds of the block ! 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 else ! allocate write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Cannot allocate meta block!" status = 1 end if ! allocate !------------------------------------------------------------------------------- ! end subroutine allocate_metablock ! !=============================================================================== ! ! 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; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine deallocate_metablock(pmeta, status) use iso_fortran_env, only : error_unit implicit none type(block_meta), pointer, intent(inout) :: pmeta integer , intent(out) :: status integer :: n, i, j #if NDIMS == 3 integer :: k #endif /* NDIMS == 3 */ ! local parameters ! character(len=*), parameter :: loc = 'BLOCKS::deallocate_metablock()' !------------------------------------------------------------------------------- ! ! check if the pointer is actually associated with any block ! if (associated(pmeta)) then ! 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 ! do i = 1, nchildren nullify(pmeta%child(i)%ptr) end do ! nullify fields pointing to face, edge, and corner neighbors ! #if NDIMS == 2 do i = 1, nsides do j = 1, nsides do n = 1, ndims nullify(pmeta%edges(i,j,n)%ptr) end do ! ndims nullify(pmeta%corners(i,j)%ptr) end do ! nsides end do ! nsides #endif /* NDIMS == 2 */ #if NDIMS == 3 do i = 1, nsides do j = 1, nsides do k = 1, nsides do n = 1, ndims nullify(pmeta%faces(i,j,k,n)%ptr) nullify(pmeta%edges(i,j,k,n)%ptr) end do ! ndims nullify(pmeta%corners(i,j,k)%ptr) end do ! nsides end do ! nsides end do ! nsides #endif /* NDIMS == 3 */ ! if there is a data block is associated, remove it ! if (associated(pmeta%data)) call remove_datablock(pmeta%data, status) ! nullify the field pointing to the associated data block ! nullify(pmeta%data) ! release the memory occupied by the block ! deallocate(pmeta, stat = status) ! nullify the pointer to the deallocated meta block ! nullify(pmeta) else ! the argument contains a null pointer, so print an error ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Null pointer argument to meta block!" status = 1 end if !------------------------------------------------------------------------------- ! end subroutine deallocate_metablock ! !=============================================================================== ! ! 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; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine allocate_datablock(pdata, status) use iso_fortran_env, only : error_unit implicit none type(block_data), pointer, intent(out) :: pdata integer , intent(out) :: status character(len=*), parameter :: loc = 'BLOCKS::allocate_datablock()' !------------------------------------------------------------------------------- ! ! allocate the block structure ! allocate(pdata, stat = status) ! check if allocation succeeded ! if (status == 0) then ! 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 and primitive variables, and fluxes ! allocate(pdata%uu(nvars,nx,ny,nz,nregs), pdata%du(nvars,nx,ny,nz), & pdata%q(nvars,nx,ny,nz), & pdata%fx(nflux,ny,nz,2), pdata%fy(nflux,nx,nz,2), & #if NDIMS == 3 pdata%fz(nflux,nx,ny,2), & #endif /* NDIMS == 3 */ stat = status) ! check if allocation succeeded ! if (status == 0) then ! reset arrays ! pdata%uu = 0.0d+00 pdata%du = 0.0d+00 pdata%q = 0.0d+00 pdata%fx = 0.0d+00 pdata%fy = 0.0d+00 #if NDIMS == 3 pdata%fz = 0.0d+00 #endif /* NDIMS == 3 */ ! initiate the conserved variable pointer ! pdata%u => pdata%uu(:,:,:,:,1) else ! allocation write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Cannot allocate variables in data block!" status = 1 end if ! allocation else ! allocation write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Cannot allocate data block!" status = 1 end if ! allocation !------------------------------------------------------------------------------- ! 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; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine deallocate_datablock(pdata, status) use iso_fortran_env, only : error_unit implicit none type(block_data), pointer, intent(inout) :: pdata integer , intent(out) :: status character(len=*), parameter :: loc = 'BLOCKS::deallocate_datablock()' !------------------------------------------------------------------------------- ! ! check if the pointer is actually associated with any block ! if (associated(pdata)) then ! 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%uu)) deallocate(pdata%uu, stat = status) ! deallocate the increment of conserved variables ! if (allocated(pdata%du)) deallocate(pdata%du, stat = status) ! deallocate primitive variables ! if (allocated(pdata%q )) deallocate(pdata%q , stat = status) ! deallocate numerical fluxes ! if (allocated(pdata%fx)) deallocate(pdata%fx, stat = status) if (allocated(pdata%fy)) deallocate(pdata%fy, stat = status) if (allocated(pdata%fz)) deallocate(pdata%fz, stat = status) ! release the memory occupied by the block ! deallocate(pdata, stat = status) ! nullify the pointer to the deallocated meta block ! nullify(pdata) else ! the argument contains a null pointer, so print an error ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Null pointer argument to data block!" status = 1 end if ! pdata associated with a data block !------------------------------------------------------------------------------- ! end subroutine deallocate_datablock ! !=============================================================================== ! ! 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) implicit none type(block_meta), pointer, intent(inout) :: pmeta type(block_data), pointer, intent(inout) :: pdata !------------------------------------------------------------------------------- ! 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) implicit none type(block_meta), pointer, intent(inout) :: pmeta type(block_data), pointer, intent(inout) :: pdata !------------------------------------------------------------------------------- ! nullify(pmeta%data) nullify(pdata%meta) !------------------------------------------------------------------------------- ! end subroutine unlink_blocks ! !=============================================================================== ! ! subroutine REFINE_BLOCK: ! ----------------------- ! ! Subroutine creates children of the current block and initializes their ! configuration, pointers and fields. ! ! Arguments: ! ! pmeta - a pointer to meta block for which children will be created; ! fdata - a flag indicating if data blocks for children should be allocated; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine refine_block(pmeta, fdata, status) use iso_fortran_env, only : error_unit implicit none type(block_meta), pointer, intent(inout) :: pmeta logical , intent(in) :: fdata integer , intent(out) :: status type(block_meta), pointer :: pnext, pneigh, pchild type(block_data), pointer :: pdata logical, save :: first = .true. !$omp threadprivate(first) integer :: p, q, np integer :: i, ip, ir integer :: j, jp, jr integer :: k, kp, kr real(kind=8) :: xln, yln, zln, xmn, xmx, ymn, ymx, zmn, zmx integer, dimension(0:79,nchildren), save :: order integer, dimension(0:79,nchildren), save :: config character(len=*), parameter :: loc = 'BLOCKS::refine_block()' !------------------------------------------------------------------------------- ! status = 0 #if NDIMS == 2 kp = 1 kr = 1 #endif /* NDIMS == 2 */ ! prepare some arrays ! if (first) then ! prepare order array ! do p = 1, nchildren order ( :,p) = p end do #if NDIMS == 2 order ( 0,:) = (/ 1, 2, 3, 4 /) order (12,:) = (/ 1, 3, 4, 2 /) order (13,:) = (/ 1, 2, 4, 3 /) order (42,:) = (/ 4, 3, 1, 2 /) order (43,:) = (/ 4, 2, 1, 3 /) #endif /* NDIMS == 2 */ #if NDIMS == 3 order ( 0,:) = (/ 1, 2, 3, 4, 5, 6, 7, 8 /) order (12,:) = (/ 1, 3, 7, 5, 6, 8, 4, 2 /) order (13,:) = (/ 1, 5, 6, 2, 4, 8, 7, 3 /) order (15,:) = (/ 1, 2, 4, 3, 7, 8, 6, 5 /) order (42,:) = (/ 4, 8, 7, 3, 1, 5, 6, 2 /) order (43,:) = (/ 4, 2, 6, 8, 7, 5, 1, 3 /) order (48,:) = (/ 4, 3, 1, 2, 6, 5, 7, 8 /) order (62,:) = (/ 6, 5, 7, 8, 4, 3, 1, 2 /) order (65,:) = (/ 6, 8, 4, 2, 1, 3, 7, 5 /) order (68,:) = (/ 6, 2, 1, 5, 7, 3, 4 , 8 /) order (73,:) = (/ 7, 8, 6, 5, 1, 2, 4, 3 /) order (75,:) = (/ 7, 3, 4, 8, 6, 2, 1, 5 /) order (78,:) = (/ 7, 5, 1, 3, 4, 2, 6, 8 /) #endif /* NDIMS == 3 */ ! prepare config array ! config( :,:) = 0 #if NDIMS == 2 config( 0,:) = (/ 0, 0, 0, 0 /) config(12,:) = (/ 13, 12, 12, 42 /) config(13,:) = (/ 12, 13, 13, 43 /) config(42,:) = (/ 43, 42, 42, 12 /) config(43,:) = (/ 42, 43, 43, 13 /) #endif /* NDIMS == 2 */ #if NDIMS == 3 config( 0,:) = (/ 0, 0, 0, 0, 0, 0, 0, 0 /) config(12,:) = (/ 13, 15, 15, 78, 78, 62, 62, 42 /) config(13,:) = (/ 15, 12, 12, 68, 68, 43, 43, 73 /) config(15,:) = (/ 12, 13, 13, 48, 48, 75, 75, 65 /) config(42,:) = (/ 48, 43, 43, 75, 75, 12, 12, 62 /) config(43,:) = (/ 42, 48, 48, 65, 65, 73, 73, 13 /) config(48,:) = (/ 43, 42, 42, 15, 15, 68, 68, 78 /) config(62,:) = (/ 65, 68, 68, 73, 73, 42, 42, 12 /) config(65,:) = (/ 68, 62, 62, 43, 43, 15, 15, 75 /) config(68,:) = (/ 62, 65, 65, 13, 13, 78, 78, 48 /) config(73,:) = (/ 78, 75, 75, 62, 62, 13, 13, 43 /) config(75,:) = (/ 73, 78, 78, 42, 42, 65, 65, 15 /) config(78,:) = (/ 75, 73, 73, 12, 12, 48, 48, 68 /) #endif /* NDIMS == 3 */ ! reset the first execution flag ! first = .false. end if ! check if pointer is associated ! if (associated(pmeta)) then ! store the pointer to the next block on the list ! pnext => pmeta%next !! PREPARE CHILD CONFIGURATION PARAMETERS !! ! set corresponding configuration of the new blocks ! q = pmeta%conf ! calculate sizes of the child blocks ! xln = 0.5d+00 * (pmeta%xmax - pmeta%xmin) yln = 0.5d+00 * (pmeta%ymax - pmeta%ymin) #if NDIMS == 3 zln = 0.5d+00 * (pmeta%zmax - pmeta%zmin) #else /* NDIMS == 3 */ zln = (pmeta%zmax - pmeta%zmin) #endif /* NDIMS == 3 */ !! ALLOCATE CHILDREN AND APPEND THEM TO THE META LIST !! ! iterate over the number of children in the reverse configuration order, i.e. ! the allocated blocks are inserted after the parent block following ! the reversed Hilbert curve ! p = nchildren do while(p >= 1 .and. status == 0) ! insert a new meta block after pmeta and associate it with pchild ! call insert_metablock_after(pmeta, pchild, status) ! proceed, if metablock insertion was successful ! if (status == 0) then ! set the child configuration number ! call metablock_set_configuration(pchild, config(q,p)) ! associate the parent's children array element with the freshly created ! meta block ! pmeta%child(order(q,p))%ptr => pchild else write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Cannot insert new meta block for the current child!" status = 1 end if ! decrease the child index ! p = p - 1 end do ! nchildren ! check if new child blocks are ready for configuration ! if (status == 0) then ! iterate over all children ! do p = 1, nchildren ! associate a pointer with the current child ! pchild => pmeta%child(p)%ptr ! associate the parent field with pmeta ! pchild%parent => pmeta ! mark the child as the leaf ! call metablock_set_leaf(pchild) ! mark the child to be updated ! call metablock_set_update(pchild) ! set the child refinement level ! call metablock_set_level(pchild, pmeta%level + 1) ! set the child process number ! call metablock_set_process(pchild, pmeta%process) ! calculate the block position indices ! q = p - 1 i = mod(q ,2) j = mod(q / 2,2) k = mod(q / 4,2) ! calculate the block coordinates in effective resolution units ! ip = 2 * pmeta%coords(1) + i jp = 2 * pmeta%coords(2) + j #if NDIMS == 3 kp = 2 * pmeta%coords(3) + k #endif /* NDIMS == 3 */ ! calculate block bounds ! xmn = pmeta%xmin + xln * i ymn = pmeta%ymin + yln * j zmn = pmeta%zmin + zln * k xmx = xmn + xln ymx = ymn + yln zmx = zmn + zln ! set the block position ! call metablock_set_position(pchild, (/ i, j, k /)) ! set the effective resolution coordinates ! call metablock_set_coordinates(pchild, (/ ip, jp, kp /)) ! set the child block bounds ! call metablock_set_bounds(pchild, xmn, xmx, ymn, ymx, zmn, zmx) end do ! nchildren ! update neighbor pointers of the parent block ! #if NDIMS == 2 do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip ! calculate the child index ! p = 2 * (jp - 1) + ip ! associate pchild with the proper child ! pchild => pmeta%child(p)%ptr !--- update edge neighbor pointers --- ! ! update external edges ! ! along X-direction ! pneigh => pmeta%edges(ip,jp,1)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 2 * (jr - 1) + ip pchild%edges(ip,jp,1)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jp,1)%ptr => pmeta%child(q)%ptr else pchild%edges(ip,jp,1)%ptr => pneigh pchild%edges(ir,jp,1)%ptr => pneigh end if end if ! pneigh associated ! along Y-direction ! pneigh => pmeta%edges(ip,jp,2)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 2 * (jp - 1) + ir pchild%edges(ip,jp,2)%ptr => pmeta%child(q)%ptr pchild%edges(ip,jr,2)%ptr => pmeta%child(q)%ptr else pchild%edges(ip,jp,2)%ptr => pneigh pchild%edges(ip,jr,2)%ptr => pneigh end if end if ! pneigh associated ! update internal edges ! ! along X-direction ! q = 2 * (jr - 1) + ip pchild%edges(ip,jr,1)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jr,1)%ptr => pmeta%child(q)%ptr ! along Y-direction ! q = 2 * (jp - 1) + ir pchild%edges(ir,jp,2)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jr,2)%ptr => pmeta%child(q)%ptr !--- update corner neighbor pointers --- ! ! calculate the index of opposite child ! q = 2 * (jr - 1) + ir ! update corner located at the parent's one ! pneigh => pmeta%corners(ip,jp)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ip,jp)%ptr => pmeta%child(q)%ptr else pchild%corners(ip,jp)%ptr => pneigh end if end if ! pneigh associated ! update corner touching another child ! pchild%corners(ir,jr)%ptr => pmeta%child(q)%ptr ! update corners laying on parent's edges ! ! along X-direction ! pneigh => pmeta%edges(ir,jp,1)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ir,jp)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) & pchild%corners(ir,jp)%ptr => pneigh end if end if ! pneigh associated ! along Y-direction ! pneigh => pmeta%edges(ip,jr,2)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ip,jr)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) & pchild%corners(ip,jr)%ptr => pneigh end if end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides #endif /* NDIMS == 2 */ #if NDIMS == 3 do kp = 1, nsides kr = 3 - kp do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip ! calculate the child index ! p = 4 * (kp - 1) + 2 * (jp - 1) + ip ! associate pchild with the proper child ! pchild => pmeta%child(p)%ptr !--- update face neighbor pointers --- ! ! prepare the index of neighbor child for X-faces ! q = 4 * (kp - 1) + 2 * (jp - 1) + ir ! set the internal side neighbor pointer ! do k = 1, nsides do j = 1, nsides pchild%faces(ir,j,k,1)%ptr => pmeta%child(q)%ptr end do ! j = 1, nsides end do ! k = 1, nsides ! associate pneigh with the X-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,1)%ptr ! set the external side neighbor pointer ! if (associated(pneigh)) then if (pneigh%id == pmeta%id) then do k = 1, nsides do j = 1, nsides pchild%faces(ip,j,k,1)%ptr => pmeta%child(q)%ptr end do ! j = 1, nsides end do ! k = 1, nsides else do k = 1, nsides do j = 1, nsides pchild%faces(ip,j,k,1)%ptr => pneigh end do ! j = 1, nsides end do ! k = 1, nsides end if end if ! pneigh associated ! prepare the index of neighbor child for Y-faces ! q = 4 * (kp - 1) + 2 * (jr - 1) + ip ! set the internal side neighbor pointer ! do k = 1, nsides do i = 1, nsides pchild%faces(i,jr,k,2)%ptr => pmeta%child(q)%ptr end do ! i = 1, nsides end do ! k = 1, nsides ! associate pneigh with the Y-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,2)%ptr ! set the external side neighbor pointer ! if (associated(pneigh)) then if (pneigh%id == pmeta%id) then do k = 1, nsides do i = 1, nsides pchild%faces(i,jp,k,2)%ptr => pmeta%child(q)%ptr end do ! i = 1, nsides end do ! k = 1, nsides else do k = 1, nsides do i = 1, nsides pchild%faces(i,jp,k,2)%ptr => pneigh end do ! i = 1, nsides end do ! k = 1, nsides end if end if ! pneigh associated ! prepare the index of neighbor child for Z-faces ! q = 4 * (kr - 1) + 2 * (jp - 1) + ip ! set the internal side neighbor pointer ! do j = 1, nsides do i = 1, nsides pchild%faces(i,j,kr,3)%ptr => pmeta%child(q)%ptr end do ! i = 1, nsides end do ! j = 1, nsides ! associate pneigh with the Z-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,3)%ptr ! set the external side neighbor pointer ! if (associated(pneigh)) then if (pneigh%id == pmeta%id) then do j = 1, nsides do i = 1, nsides pchild%faces(i,j,kp,3)%ptr => pmeta%child(q)%ptr end do ! i = 1, nsides end do ! j = 1, nsides else do j = 1, nsides do i = 1, nsides pchild%faces(i,j,kp,3)%ptr => pneigh end do ! i = 1, nsides end do ! j = 1, nsides end if end if ! pneigh associated !--- update edge neighbor pointers --- ! ! process child edges which lay on the parent's edges ! ! along X direction ! pneigh => pmeta%edges(ip,jp,kp,1)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kr - 1) + 2 * (jr - 1) + ip pchild%edges(ip,jp,kp,1)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jp,kp,1)%ptr => pmeta%child(q)%ptr else pchild%edges(ip,jp,kp,1)%ptr => pneigh pchild%edges(ir,jp,kp,1)%ptr => pneigh end if end if ! pneigh associated ! along Y direction ! pneigh => pmeta%edges(ip,jp,kp,2)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kr - 1) + 2 * (jp - 1) + ir pchild%edges(ip,jp,kp,2)%ptr => pmeta%child(q)%ptr pchild%edges(ip,jr,kp,2)%ptr => pmeta%child(q)%ptr else pchild%edges(ip,jp,kp,2)%ptr => pneigh pchild%edges(ip,jr,kp,2)%ptr => pneigh end if end if ! pneigh associated ! along Z direction ! pneigh => pmeta%edges(ip,jp,kp,3)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kp - 1) + 2 * (jr - 1) + ir pchild%edges(ip,jp,kp,3)%ptr => pmeta%child(q)%ptr pchild%edges(ip,jp,kr,3)%ptr => pmeta%child(q)%ptr else pchild%edges(ip,jp,kp,3)%ptr => pneigh pchild%edges(ip,jp,kr,3)%ptr => pneigh end if end if ! pneigh associated ! process child edges which are neighbors with other children ! ! along X direction ! q = 4 * (kr - 1) + 2 * (jr - 1) + ip pchild%edges(ip,jr,kr,1)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jr,kr,1)%ptr => pmeta%child(q)%ptr ! along Y direction ! q = 4 * (kr - 1) + 2 * (jp - 1) + ir pchild%edges(ir,jp,kr,2)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jr,kr,2)%ptr => pmeta%child(q)%ptr ! along Z direction ! q = 4 * (kp - 1) + 2 * (jr - 1) + ir pchild%edges(ir,jr,kp,3)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jr,kr,3)%ptr => pmeta%child(q)%ptr ! process child edges on the parent's X-face ! ! along Z-edge ! pneigh => pmeta%faces(ip,jr,kp,1)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kp - 1) + 2 * (jr - 1) + ir pchild%edges(ip,jr,kp,3)%ptr => pmeta%child(q)%ptr pchild%edges(ip,jr,kr,3)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) then pchild%edges(ip,jr,kp,3)%ptr => pneigh pchild%edges(ip,jr,kr,3)%ptr => pneigh end if end if end if ! pneigh associated ! along Y-edge ! pneigh => pmeta%faces(ip,jp,kr,1)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kr - 1) + 2 * (jp - 1) + ir pchild%edges(ip,jp,kr,2)%ptr => pmeta%child(q)%ptr pchild%edges(ip,jr,kr,2)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) then pchild%edges(ip,jp,kr,2)%ptr => pneigh pchild%edges(ip,jr,kr,2)%ptr => pneigh end if end if end if ! pneigh associated ! process child edges on the parent's Y-face ! ! along Z-edge ! pneigh => pmeta%faces(ir,jp,kp,2)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kp - 1) + 2 * (jr - 1) + ir pchild%edges(ir,jp,kp,3)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jp,kr,3)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) then pchild%edges(ir,jp,kp,3)%ptr => pneigh pchild%edges(ir,jp,kr,3)%ptr => pneigh end if end if end if ! pneigh associated ! along X-edge ! pneigh => pmeta%faces(ip,jp,kr,2)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kr - 1) + 2 * (jr - 1) + ip pchild%edges(ip,jp,kr,1)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jp,kr,1)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) then pchild%edges(ip,jp,kr,1)%ptr => pneigh pchild%edges(ir,jp,kr,1)%ptr => pneigh end if end if end if ! pneigh associated ! process child edges on the parent's Z-face ! ! along Y-edge ! pneigh => pmeta%faces(ir,jp,kp,3)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kr - 1) + 2 * (jp - 1) + ir pchild%edges(ir,jp,kp,2)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jr,kp,2)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) then pchild%edges(ir,jp,kp,2)%ptr => pneigh pchild%edges(ir,jr,kp,2)%ptr => pneigh end if end if end if ! pneigh associated ! along X-edge ! pneigh => pmeta%faces(ip,jr,kp,3)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then q = 4 * (kr - 1) + 2 * (jr - 1) + ip pchild%edges(ip,jr,kp,1)%ptr => pmeta%child(q)%ptr pchild%edges(ir,jr,kp,1)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) then pchild%edges(ip,jr,kp,1)%ptr => pneigh pchild%edges(ir,jr,kp,1)%ptr => pneigh end if end if end if ! pneigh associated !--- update corner neighbor pointers --- ! ! calculate the index of the neighbor child ! q = 4 * (kr - 1) + 2 * (jr - 1) + ir ! process child corner which overlaps with the parent's one ! pneigh => pmeta%corners(ip,jp,kp)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ip,jp,kp)%ptr => pmeta%child(q)%ptr else pchild%corners(ip,jp,kp)%ptr => pmeta%corners(ip,jp,kp)%ptr end if end if ! pneigh associated ! process child corner which points to another child ! pchild%corners(ir,jr,kr)%ptr => pmeta%child(q)%ptr ! process child corners which lay on parent's edges ! ! along X direction ! pneigh => pmeta%edges(ir,jp,kp,1)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ir,jp,kp)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) & pchild%corners(ir,jp,kp)%ptr => pneigh end if end if ! pneigh associated ! along Y direction ! pneigh => pmeta%edges(ip,jr,kp,2)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ip,jr,kp)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) & pchild%corners(ip,jr,kp)%ptr => pneigh end if end if ! pneigh associated ! along Z-direction ! pneigh => pmeta%edges(ip,jp,kr,3)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ip,jp,kr)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) & pchild%corners(ip,jp,kr)%ptr => pneigh end if end if ! pneigh associated ! process child corners which lay on parent's faces ! ! on X-face ! pneigh => pmeta%faces(ip,jr,kr,1)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ip,jr,kr)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) & pchild%corners(ip,jr,kr)%ptr => pneigh end if end if ! pneigh associated ! on Y-face ! pneigh => pmeta%faces(ir,jp,kr,2)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ir,jp,kr)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) & pchild%corners(ir,jp,kr)%ptr => pneigh end if end if ! pneigh associated ! on Z-face ! pneigh => pmeta%faces(ir,jr,kp,3)%ptr if (associated(pneigh)) then if (pneigh%id == pmeta%id) then pchild%corners(ir,jr,kp)%ptr => pmeta%child(q)%ptr else if (pneigh%level > pmeta%level) & pchild%corners(ir,jr,kp)%ptr => pneigh end if end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides end do ! kp = 1, nsides #endif /* NDIMS == 3 */ ! update neighbor pointers of the neighbor blocks ! #if NDIMS == 2 do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip ! calculate the child index ! p = 2 * (jp - 1) + ip ! associate pchild with the proper child ! pchild => pmeta%child(p)%ptr !--- update neighbor's edge pointers --- ! ! along X-direction ! pneigh => pchild%edges(ip,jp,1)%ptr if (associated(pneigh)) then pneigh%edges(ip,jr,1)%ptr => pchild if (pneigh%level > pmeta%level) & pneigh%edges(ir,jr,1)%ptr => pchild end if ! pneigh associated ! along Y-direction ! pneigh => pchild%edges(ip,jp,2)%ptr if (associated(pneigh)) then pneigh%edges(ir,jp,2)%ptr => pchild if (pneigh%level > pmeta%level) & pneigh%edges(ir,jr,2)%ptr => pchild end if ! pneigh associated !--- update neighbor's corner pointers --- ! ! neighbor corner located at the parent's one ! pneigh => pmeta%corners(ip,jp)%ptr if (associated(pneigh)) pneigh%corners(ir,jr)%ptr => pchild ! neighbor corners laying on the parent's edges ! ! along X-direction ! pneigh => pmeta%edges(ir,jp,1)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & pneigh%corners(ip,jr)%ptr => pchild end if ! pneigh associated ! along Y-direction ! pneigh => pmeta%edges(ip,jr,2)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & pneigh%corners(ir,jp)%ptr => pchild end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides #endif /* NDIMS == 2 */ #if NDIMS == 3 ! update neighbor's face pointers (only in 3D) ! do kp = 1, nsides kr = 3 - kp do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip ! calculate the child index ! p = 4 * (kp - 1) + 2 * (jp - 1) + ip ! associate pchild with the proper child ! pchild => pmeta%child(p)%ptr !--- update neighbor's face pointers --- ! ! assign pneigh to the X-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,1)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then if (pneigh%level > pmeta%level) then do k = 1, nsides do j = 1, nsides pneigh%faces(ir,j,k,1)%ptr => pchild end do end do else pneigh%faces(ir,jp,kp,1)%ptr => pchild end if end if ! pneigh associated ! assign pneigh to the Y-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,2)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then if (pneigh%level > pmeta%level) then do k = 1, nsides do i = 1, nsides pneigh%faces(i,jr,k,2)%ptr => pchild end do end do else pneigh%faces(ip,jr,kp,2)%ptr => pchild end if end if ! pneigh associated ! assign pneigh to the Z-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,3)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then if (pneigh%level > pmeta%level) then do j = 1, nsides do i = 1, nsides pneigh%faces(i,j,kr,3)%ptr => pchild end do end do else pneigh%faces(ip,jp,kr,3)%ptr => pchild end if end if ! pneigh associated !--- update neighbor's edge pointers --- ! ! along X direction ! pneigh => pmeta%edges(ip,jp,kp,1)%ptr if (associated(pneigh)) then pneigh%edges(ip,jr,kr,1)%ptr => pchild if (pneigh%level > pmeta%level) & pneigh%edges(ir,jr,kr,1)%ptr => pchild end if ! pneigh associated ! along Y direction ! pneigh => pmeta%edges(ip,jp,kp,2)%ptr if (associated(pneigh)) then pneigh%edges(ir,jp,kr,2)%ptr => pchild if (pneigh%level > pmeta%level) & pneigh%edges(ir,jr,kr,2)%ptr => pchild end if ! pneigh associated ! along Z direction ! pneigh => pmeta%edges(ip,jp,kp,3)%ptr if (associated(pneigh)) then pneigh%edges(ir,jr,kp,3)%ptr => pchild if (pneigh%level > pmeta%level) & pneigh%edges(ir,jr,kr,3)%ptr => pchild end if ! pneigh associated ! process child edges on the parent's X face ! ! Z-direction ! pneigh => pmeta%faces(ip,jr,kp,1)%ptr if (associated(pneigh)) then if (pneigh%level == pchild%level) then pneigh%edges(ir,jp,kp,3)%ptr => pchild pneigh%edges(ir,jp,kr,3)%ptr => pchild end if end if ! pneigh associated ! Y-direction ! pneigh => pmeta%faces(ip,jp,kr,1)%ptr if (associated(pneigh)) then if (pneigh%level == pchild%level) then pneigh%edges(ir,jp,kp,2)%ptr => pchild pneigh%edges(ir,jr,kp,2)%ptr => pchild end if end if ! pneigh associated ! process child edges on the parent's Y face ! ! Z-direction ! pneigh => pmeta%faces(ir,jp,kp,2)%ptr if (associated(pneigh)) then if (pneigh%level == pchild%level) then pneigh%edges(ip,jr,kp,3)%ptr => pchild pneigh%edges(ip,jr,kr,3)%ptr => pchild end if end if ! pneigh associated ! X-direction ! pneigh => pmeta%faces(ip,jp,kr,2)%ptr if (associated(pneigh)) then if (pneigh%level == pchild%level) then pneigh%edges(ip,jr,kp,1)%ptr => pchild pneigh%edges(ir,jr,kp,1)%ptr => pchild end if end if ! pneigh associated ! process child edges on the parent's Z face ! ! Y-direction ! pneigh => pmeta%faces(ir,jp,kp,3)%ptr if (associated(pneigh)) then if (pneigh%level == pchild%level) then pneigh%edges(ip,jp,kr,2)%ptr => pchild pneigh%edges(ip,jr,kr,2)%ptr => pchild end if end if ! pneigh associated ! X-direction ! pneigh => pmeta%faces(ip,jr,kp,3)%ptr if (associated(pneigh)) then if (pneigh%level == pchild%level) then pneigh%edges(ip,jp,kr,1)%ptr => pchild pneigh%edges(ir,jp,kr,1)%ptr => pchild end if end if ! pneigh associated !--- update neighbor's corner pointers --- ! ! calculate the index of the opposite child ! q = 4 * (kr - 1) + 2 * (jr - 1) + ir ! update neighbor's corner which overlaps with the parent's one ! pneigh => pmeta%corners(ip,jp,kp)%ptr if (associated(pneigh)) pneigh%corners(ir,jr,kr)%ptr => pchild ! process neighbot's corners which lay on parent's edges ! ! X-edge ! pneigh => pmeta%edges(ir,jp,kp,1)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & pneigh%corners(ip,jr,kr)%ptr => pchild end if ! pneigh associated ! Y-edge ! pneigh => pmeta%edges(ip,jr,kp,2)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & pneigh%corners(ir,jp,kr)%ptr => pchild end if ! pneigh associated ! Z-edge ! pneigh => pmeta%edges(ip,jp,kr,3)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & pneigh%corners(ir,jr,kp)%ptr => pchild end if ! pneigh associated ! process child corners which lay on parent's faces ! ! X-face ! pneigh => pmeta%faces(ip,jr,kr,1)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & pneigh%corners(ir,jp,kp)%ptr => pchild end if ! pneigh associated ! Y-face ! pneigh => pmeta%faces(ir,jp,kr,2)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & pneigh%corners(ip,jr,kp)%ptr => pchild end if ! pneigh associated ! Z-face ! pneigh => pmeta%faces(ir,jr,kp,3)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & pneigh%corners(ip,jp,kr)%ptr => pchild end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides end do ! kp = 1, nsides #endif /* NDIMS == 3 */ !! ASSOCIATE DATA BLOCKS IF NECESSARY !! ! allocate data blocks if requested ! if (fdata) then ! iterate over all children ! p = 1 do while(p <= nchildren .and. status == 0) ! assign a pointer to the current child ! pchild => pmeta%child(p)%ptr ! allocate new data block and append it to the data block list ! call append_datablock(pdata, status) ! check if data block appended successfully ! if (status == 0) then ! associate the new data block with the current child ! call link_blocks(pchild, pdata) else write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Cannot append new data block for the current child!" status = 1 end if ! increase child index ! p = p + 1 end do ! nchildren end if ! allocate data blocks for children !! RESET PARENT'S FIELDS !! if (status == 0) then ! derefine all neighbors ! #if NDIMS == 3 do kp = 1, nsides #endif /* NDIMS == 3 */ do jp = 1, nsides do ip = 1, nsides do np = 1, ndims #if NDIMS == 2 if (associated(pmeta%edges(ip,jp,np)%ptr)) & nullify(pmeta%edges(ip,jp,np)%ptr) #endif /* NDIMS == 2 */ #if NDIMS == 3 if (associated(pmeta%faces(ip,jp,kp,np)%ptr)) & nullify(pmeta%faces(ip,jp,kp,np)%ptr) if (associated(pmeta%edges(ip,jp,kp,np)%ptr)) & nullify(pmeta%edges(ip,jp,kp,np)%ptr) #endif /* NDIMS == 3 */ end do ! np = 1, ndims #if NDIMS == 2 if (associated(pmeta%corners(ip,jp)%ptr)) & nullify(pmeta%corners(ip,jp)%ptr) #endif /* NDIMS == 2 */ #if NDIMS == 3 if (associated(pmeta%corners(ip,jp,kp)%ptr)) & nullify(pmeta%corners(ip,jp,kp)%ptr) #endif /* NDIMS == 3 */ end do ! ip = 1, nsides end do ! jp = 1, nsides #if NDIMS == 3 end do ! kp = 1, nsides #endif /* NDIMS == 3 */ ! unset the block leaf flag ! call metablock_unset_leaf(pmeta) ! reset the refinement flag ! call metablock_set_refinement(pmeta, 0) ! restore the pointer to the current block ! if (associated(pnext)) then pmeta => pnext%prev else pmeta => last_meta end if end if ! status end if ! children ready else ! pmeta is not associated write(error_unit,"('[',a,']: ',a)") trim(loc) & , "No meta block associated with the pointer argument!" status = 1 end if !------------------------------------------------------------------------------- ! end subroutine refine_block ! !=============================================================================== ! ! subroutine DEREFINE_BLOCK: ! ------------------------- ! ! Subroutine derefines the current block by distrying all its children and ! restoring the block configuration, pointers and fields. ! ! Arguments: ! ! pmeta - a pointer to derefined meta block; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine derefine_block(pmeta, status) implicit none type(block_meta), pointer, intent(inout) :: pmeta integer , intent(out) :: status type(block_meta), pointer :: pchild, pneigh integer :: p, q integer :: ip, ir integer :: jp, jr #if NDIMS == 3 integer :: i, j, k integer :: kp, kr #endif /* NDIMS == 3 */ !------------------------------------------------------------------------------- ! status = 0 ! update neighbor pointers of the parent block ! #if NDIMS == 2 do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip ! calculate the child index ! p = 2 * (jp - 1) + ip ! associate pchild with the proper child ! pchild => pmeta%child(p)%ptr !--- update edge neighbor pointers --- ! ! along X-direction ! pneigh => pchild%edges(ip,jp,1)%ptr if (associated(pneigh)) then q = 2 * (jr - 1) + ip if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%edges(ip,jp,1)%ptr => pmeta else pmeta%edges(ip,jp,1)%ptr => pneigh end if end if ! pneigh associated ! along Y-direction ! pneigh => pchild%edges(ip,jp,2)%ptr if (associated(pneigh)) then q = 2 * (jp - 1) + ir if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%edges(ip,jp,2)%ptr => pmeta else pmeta%edges(ip,jp,2)%ptr => pneigh end if end if ! pneigh associated !--- update corner neighbor pointers --- ! pneigh => pchild%corners(ip,jp)%ptr if (associated(pneigh)) then q = 2 * (jr - 1) + ir if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%corners(ip,jp)%ptr => pmeta else pmeta%corners(ip,jp)%ptr => pneigh end if end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides #endif /* NDIMS == 2 */ #if NDIMS == 3 do kp = 1, nsides kr = 3 - kp do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip ! calculate the child index ! p = 4 * (kp - 1) + 2 * (jp - 1) + ip ! associate pchild with the proper child ! pchild => pmeta%child(p)%ptr !--- update face neighbor pointers --- ! ! assign pneigh to the X-face neighbor ! pneigh => pchild%faces(ip,jp,kp,1)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then q = 4 * (kp - 1) + 2 * (jp - 1) + ir if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%faces(ip,jp,kp,1)%ptr => pmeta else pmeta%faces(ip,jp,kp,1)%ptr => pneigh end if end if ! assign pneigh to the Y-face neighbor ! pneigh => pchild%faces(ip,jp,kp,2)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then q = 4 * (kp - 1) + 2 * (jr - 1) + ip if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%faces(ip,jp,kp,2)%ptr => pmeta else pmeta%faces(ip,jp,kp,2)%ptr => pneigh end if end if ! assign pneigh to the Z-face neighbor ! pneigh => pchild%faces(ip,jp,kp,3)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then q = 4 * (kr - 1) + 2 * (jp - 1) + ip if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%faces(ip,jp,kp,3)%ptr => pmeta else pmeta%faces(ip,jp,kp,3)%ptr => pneigh end if end if !--- update edge neighbor pointers --- ! ! associate pneigh with the X edge neighbor ! pneigh => pchild%edges(ip,jp,kp,1)%ptr ! process edge along X-direction if pneigh associated ! if (associated(pneigh)) then q = 4 * (kr - 1) + 2 * (jr - 1) + ip if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%edges(ip,jp,kp,1)%ptr => pmeta else pmeta%edges(ip,jp,kp,1)%ptr => pneigh end if end if ! pneigh associated ! associate pneigh with the Y edge neighbor ! pneigh => pchild%edges(ip,jp,kp,2)%ptr ! process edge along Y-direction if pneigh associated ! if (associated(pneigh)) then q = 4 * (kr - 1) + 2 * (jp - 1) + ir if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%edges(ip,jp,kp,2)%ptr => pmeta else pmeta%edges(ip,jp,kp,2)%ptr => pneigh end if end if ! pneigh associated ! associate pneigh with the Z edge neighbor ! pneigh => pchild%edges(ip,jp,kp,3)%ptr ! process edge along Y-direction if pneigh associated ! if (associated(pneigh)) then q = 4 * (kp - 1) + 2 * (jr - 1) + ir if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%edges(ip,jp,kp,3)%ptr => pmeta else pmeta%edges(ip,jp,kp,3)%ptr => pneigh end if end if ! pneigh associated !--- update corner neighbor pointers --- ! ! associate pneigh with the corner neighbor ! pneigh => pchild%corners(ip,jp,kp)%ptr ! update the corner neighbor pointer ! if (associated(pneigh)) then ! calculate the index of the opposite child ! q = 4 * (kr - 1) + 2 * (jr - 1) + ir if (pneigh%id == pmeta%child(q)%ptr%id) then pmeta%corners(ip,jp,kp)%ptr => pmeta else pmeta%corners(ip,jp,kp)%ptr => pneigh end if end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides end do ! kp = 1, nsides #endif /* NDIMS == 3 */ ! update neighbor pointers of the neighbor blocks ! #if NDIMS == 2 do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip !--- update neighbor's edge pointers --- ! ! along X-direction ! pneigh => pmeta%edges(ip,jp,1)%ptr if (associated(pneigh)) then pneigh%edges(ip,jr,1)%ptr => pmeta if (pneigh%level > pmeta%level) pneigh%edges(ir,jr,1)%ptr => pmeta end if ! pneigh associated ! along Y-direction ! pneigh => pmeta%edges(ip,jp,2)%ptr if (associated(pneigh)) then pneigh%edges(ir,jp,2)%ptr => pmeta if (pneigh%level > pmeta%level) pneigh%edges(ir,jr,2)%ptr => pmeta end if !--- update neighbor's corner pointers --- ! ! neighbor corner linked to the parent's corner ! pneigh => pmeta%corners(ip,jp)%ptr if (associated(pneigh)) pneigh%corners(ir,jr)%ptr => pmeta ! nullify neighbor corners pointing to parent's edges ! ! along X-direction ! pneigh => pmeta%edges(ir,jp,1)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) nullify(pneigh%corners(ip,jr)%ptr) end if ! pneigh associated ! along Y-direction ! pneigh => pmeta%edges(ip,jr,2)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) nullify(pneigh%corners(ir,jp)%ptr) end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides #endif /* NDIMS == 2 */ #if NDIMS == 3 do kp = 1, nsides kr = 3 - kp do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip !--- update neighbor's face pointers --- ! ! assign pneigh to the X-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,1)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then if (pneigh%level > pmeta%level) then do k = 1, nsides do j = 1, nsides pneigh%faces(ir,j,k,1)%ptr => pmeta end do end do else pneigh%faces(ir,jp,kp,1)%ptr => pmeta end if end if ! pneigh associated ! assign pneigh to the Y-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,2)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then if (pneigh%level > pmeta%level) then do k = 1, nsides do i = 1, nsides pneigh%faces(i,jr,k,2)%ptr => pmeta end do end do else pneigh%faces(ip,jr,kp,2)%ptr => pmeta end if end if ! pneigh associated ! assign pneigh to the Z-face neighbor ! pneigh => pmeta%faces(ip,jp,kp,3)%ptr ! set the corresponding neighbor face pointers ! if (associated(pneigh)) then if (pneigh%level > pmeta%level) then do j = 1, nsides do i = 1, nsides pneigh%faces(i,j,kr,3)%ptr => pmeta end do end do else pneigh%faces(ip,jp,kr,3)%ptr => pmeta end if end if ! pneigh associated !--- update neighbor's edge pointers --- ! ! process the edges in all directions which lay on the parent's edges ! ! X-edge ! pneigh => pmeta%edges(ip,jp,kp,1)%ptr if (associated(pneigh)) then pneigh%edges(ip,jr,kr,1)%ptr => pmeta if (pneigh%level > pmeta%level) & pneigh%edges(ir,jr,kr,1)%ptr => pmeta end if ! pneigh associated ! Y-edge ! pneigh => pmeta%edges(ip,jp,kp,2)%ptr if (associated(pneigh)) then pneigh%edges(ir,jp,kr,2)%ptr => pmeta if (pneigh%level > pmeta%level) & pneigh%edges(ir,jr,kr,2)%ptr => pmeta end if ! pneigh associated ! Z-edge ! pneigh => pmeta%edges(ip,jp,kp,3)%ptr if (associated(pneigh)) then pneigh%edges(ir,jr,kp,3)%ptr => pmeta if (pneigh%level > pmeta%level) & pneigh%edges(ir,jr,kr,3)%ptr => pmeta end if ! pneigh associated ! nullify neighbors edge pointers if they are on higher levels ! ! X-face ! pneigh => pmeta%faces(ip,jr,kp,1)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) then nullify(pneigh%edges(ir,jp,kp,3)%ptr) nullify(pneigh%edges(ir,jp,kr,3)%ptr) end if end if ! pneigh associated pneigh => pmeta%faces(ip,jp,kr,1)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) then nullify(pneigh%edges(ir,jp,kp,2)%ptr) nullify(pneigh%edges(ir,jr,kp,2)%ptr) end if end if ! pneigh associated ! Y-face ! pneigh => pmeta%faces(ir,jp,kp,2)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) then nullify(pneigh%edges(ip,jr,kp,3)%ptr) nullify(pneigh%edges(ip,jr,kr,3)%ptr) end if end if ! pneigh associated pneigh => pmeta%faces(ip,jp,kr,2)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) then nullify(pneigh%edges(ip,jr,kp,1)%ptr) nullify(pneigh%edges(ir,jr,kp,1)%ptr) end if end if ! pneigh associated ! Z-face ! pneigh => pmeta%faces(ir,jp,kp,3)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) then nullify(pneigh%edges(ip,jp,kr,2)%ptr) nullify(pneigh%edges(ip,jr,kr,2)%ptr) end if end if ! pneigh associated pneigh => pmeta%faces(ip,jr,kp,3)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) then nullify(pneigh%edges(ip,jp,kr,1)%ptr) nullify(pneigh%edges(ir,jp,kr,1)%ptr) end if end if ! pneigh associated !--- update neighbor's corner pointers --- ! ! associate pneigh with the corner pointer ! pneigh => pmeta%corners(ip,jp,kp)%ptr if (associated(pneigh)) pneigh%corners(ir,jr,kr)%ptr => pmeta ! nullify neighbor corners pointing to pmeta edges ! pneigh => pmeta%edges(ir,jp,kp,1)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & nullify(pneigh%corners(ip,jr,kr)%ptr) end if ! pneigh associated pneigh => pmeta%edges(ip,jr,kp,2)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & nullify(pneigh%corners(ir,jp,kr)%ptr) end if ! pneigh associated pneigh => pmeta%edges(ip,jp,kr,3)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & nullify(pneigh%corners(ir,jr,kp)%ptr) end if ! pneigh associated ! nullify neighbor corners pointing to pmeta faces ! pneigh => pmeta%faces(ip,jr,kr,1)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & nullify(pneigh%corners(ir,jp,kp)%ptr) end if ! pneigh associated pneigh => pmeta%faces(ir,jp,kr,2)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & nullify(pneigh%corners(ip,jr,kp)%ptr) end if ! pneigh associated pneigh => pmeta%faces(ir,jr,kp,3)%ptr if (associated(pneigh)) then if (pneigh%level > pmeta%level) & nullify(pneigh%corners(ip,jp,kr)%ptr) end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides end do ! kp = 1, nsides #endif /* NDIMS == 3 */ ! remove the children from the meta block list ! p = 1 do while(p <= nchildren .and. status == 0) call remove_metablock(pmeta%child(p)%ptr, status) p = p + 1 end do ! nchild ! update the parent leaf flag ! call metablock_set_leaf(pmeta) ! reset the refinement flag of the parent block ! call metablock_set_refinement(pmeta, 0) ! mark the parent to be updated ! call metablock_set_update(pmeta) !------------------------------------------------------------------------------- ! end subroutine derefine_block ! !=============================================================================== ! ! subroutine SET_LAST_ID: ! ---------------------- ! ! Subroutine sets the last identification number. This subroutine should ! be only used when the job is resumed. ! ! Arguments: ! ! id - the identification number to set; ! !=============================================================================== ! subroutine set_last_id(id) use iso_fortran_env, only : error_unit implicit none integer(kind=4), intent(in) :: id character(len=*), parameter :: loc = 'BLOCKS::set_last_id()' !------------------------------------------------------------------------------- ! ! check if the new last_id is larger than the already existing ! if (last_id > id) then write(error_unit,"('[',a,']: ',a)") trim(loc) & , "New last_id must be larger than the old one!" else ! set the last identification number ! last_id = id end if !------------------------------------------------------------------------------- ! end subroutine set_last_id ! !=============================================================================== ! ! function GET_LAST_ID: ! -------------------- ! ! Function returns the last identification number. ! ! !=============================================================================== ! function get_last_id() result(id) implicit none integer(kind=4) :: id !------------------------------------------------------------------------------- ! id = last_id !------------------------------------------------------------------------------- ! end function get_last_id ! !=============================================================================== ! ! function GET_MBLOCKS: ! -------------------- ! ! Function returns the number of meta blocks. ! !=============================================================================== ! function get_mblocks() result(nr) implicit none integer(kind=4) :: nr !------------------------------------------------------------------------------- ! nr = mblocks !------------------------------------------------------------------------------- ! end function get_mblocks ! !=============================================================================== ! ! function GET_DBLOCKS: ! -------------------- ! ! Function returns the number of data blocks. ! !=============================================================================== ! function get_dblocks() result(nr) implicit none integer(kind=4) :: nr !------------------------------------------------------------------------------- ! nr = dblocks !------------------------------------------------------------------------------- ! end function get_dblocks ! !=============================================================================== ! ! function GET_NLEAFS: ! ------------------- ! ! Function returns the number of leafs. ! !=============================================================================== ! function get_nleafs() result(nr) implicit none integer(kind=4) :: nr !------------------------------------------------------------------------------- ! nr = nleafs !------------------------------------------------------------------------------- ! end function get_nleafs ! !=============================================================================== ! ! subroutine SET_BLOCKS_UPDATE: ! ---------------------------- ! ! Subroutine sets the update flag of all meta block in the list. ! ! Arguments: ! ! flag - the flag to be set; ! !=============================================================================== ! subroutine set_blocks_update(flag) implicit none logical, intent(in) :: flag type(block_meta), pointer :: pmeta !------------------------------------------------------------------------------- ! pmeta => list_meta do while(associated(pmeta)) pmeta%update = flag pmeta => pmeta%next end do ! meta blocks !------------------------------------------------------------------------------- ! end subroutine set_blocks_update ! !=============================================================================== ! ! subroutine CHANGE_BLOCKS_PROCESS: ! -------------------------------- ! ! Subroutine changes the process of meta blocks associated with the process ! npold to npnew. ! ! Arguments: ! ! npold - the old process number; ! npnew - the new process number; ! !=============================================================================== ! subroutine change_blocks_process(npold, npnew) implicit none integer, intent(in) :: npold, npnew type(block_meta), pointer :: pmeta !------------------------------------------------------------------------------- ! if (npnew == npold) return pmeta => list_meta do while(associated(pmeta)) if (pmeta%process == npold) pmeta%process = npnew pmeta => pmeta%next end do !------------------------------------------------------------------------------- ! end subroutine change_blocks_process ! !=============================================================================== ! ! subroutine METABLOCK_SET_ID: ! --------------------------- ! ! Subroutine sets the identification number of the meta block pointed by ! the input argument. This subroutine should be used only when resuming jobs. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! id - the identification number to set; ! !=============================================================================== ! subroutine metablock_set_id(pmeta, id) implicit none type(block_meta), pointer, intent(inout) :: pmeta integer(kind=4) , intent(in) :: id !------------------------------------------------------------------------------- ! ! set the meta block %id field ! pmeta%id = id ! check if the last identification number is smaller than id, if so set ! the value of last_id to id ! if (last_id < id) last_id = id !------------------------------------------------------------------------------- ! end subroutine metablock_set_id ! !=============================================================================== ! ! subroutine METABLOCK_SET_PROCESS: ! -------------------------------- ! ! Subroutine sets the process number of the meta block pointed by ! the input argument. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! np - the process number; ! !=============================================================================== ! subroutine metablock_set_process(pmeta, np) implicit none type(block_meta), pointer, intent(inout) :: pmeta integer(kind=4) , intent(in) :: np !------------------------------------------------------------------------------- ! pmeta%process = np !------------------------------------------------------------------------------- ! end subroutine metablock_set_process ! !=============================================================================== ! ! subroutine METABLOCK_SET_LEVEL: ! ------------------------------ ! ! Subroutine sets the refinement level number of the meta block pointed ! by the input argument. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! lv - the refinement level number; ! !=============================================================================== ! subroutine metablock_set_level(pmeta, lv) implicit none type(block_meta), pointer, intent(inout) :: pmeta integer(kind=4) , intent(in) :: lv !------------------------------------------------------------------------------- ! pmeta%level = lv !------------------------------------------------------------------------------- ! end subroutine metablock_set_level ! !=============================================================================== ! ! subroutine METABLOCK_SET_CONFIGURATION: ! -------------------------------------- ! ! Subroutine sets the children block configuration number of the meta block ! pointed by the input argument. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! cf - the configuration number; ! !=============================================================================== ! subroutine metablock_set_configuration(pmeta, cf) implicit none type(block_meta), pointer, intent(inout) :: pmeta integer(kind=4) , intent(in) :: cf !------------------------------------------------------------------------------- ! pmeta%conf = cf !------------------------------------------------------------------------------- ! end subroutine metablock_set_configuration ! !=============================================================================== ! ! subroutine METABLOCK_SET_REFINEMENT: ! ----------------------------------- ! ! Subroutine sets the refinement flag of the meta block pointed by ! the input argument. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! rf - the refinement flag; ! !=============================================================================== ! subroutine metablock_set_refinement(pmeta, rf) use iso_fortran_env, only : error_unit implicit none type(block_meta), pointer, intent(inout) :: pmeta integer(kind=4) , intent(in) :: rf character(len=*), parameter :: loc = 'BLOCKS::metablock_set_refinement()' !------------------------------------------------------------------------------- ! ! check if the refinement value is correct ! if (abs(rf) > 1) then ! print error about wrong refine flag ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "The refinement value is wrong!" else ! set the block's refinement field ! pmeta%refine = rf end if !------------------------------------------------------------------------------- ! end subroutine metablock_set_refinement ! !=============================================================================== ! ! subroutine METABLOCK_SET_POSITION: ! --------------------------------- ! ! Subroutine sets the position coordinates in the parent block of ! the meta block pointed by the input argument. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! pos - the block position coordinates; ! !=============================================================================== ! subroutine metablock_set_position(pmeta, pos) implicit none type(block_meta), pointer , intent(inout) :: pmeta integer(kind=4), dimension(3), intent(in) :: pos !------------------------------------------------------------------------------- ! pmeta%pos(1:NDIMS) = pos(1:NDIMS) !------------------------------------------------------------------------------- ! end subroutine metablock_set_position ! !=============================================================================== ! ! subroutine METABLOCK_SET_COORDINATES: ! ------------------------------------ ! ! Subroutine sets the effective resolution coordinates of the meta block ! pointed by the input argument. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! coords - the effective resolution coordinates; ! !=============================================================================== ! subroutine metablock_set_coordinates(pmeta, coords) implicit none type(block_meta), pointer , intent(inout) :: pmeta integer(kind=4), dimension(3), intent(in) :: coords !------------------------------------------------------------------------------- ! pmeta%coords(1:NDIMS) = coords(1:NDIMS) !------------------------------------------------------------------------------- ! end subroutine metablock_set_coordinates ! !=============================================================================== ! ! subroutine METABLOCK_SET_BOUNDS: ! ------------------------------- ! ! Subroutine sets the physical coordinate bounds of the meta block pointed ! by the input argument. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! xmn, xmx - the coordinate bounds along X; ! ymn, ymx - the coordinate bounds along Y; ! zmn, zmx - the coordinate bounds along Z; ! !=============================================================================== ! subroutine metablock_set_bounds(pmeta, xmn, xmx, ymn, ymx, zmn, zmx) implicit none type(block_meta), pointer, intent(inout) :: pmeta real(kind=8) , intent(in) :: xmn, xmx real(kind=8) , intent(in) :: ymn, ymx real(kind=8) , intent(in) :: zmn, zmx !------------------------------------------------------------------------------- ! pmeta%xmin = xmn pmeta%xmax = xmx pmeta%ymin = ymn pmeta%ymax = ymx pmeta%zmin = zmn pmeta%zmax = zmx !------------------------------------------------------------------------------- ! end subroutine metablock_set_bounds ! !=============================================================================== ! ! subroutine METABLOCK_SET_LEAF: ! ----------------------------- ! ! Subroutine marks the meta block pointed by the input argument as the leaf. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! !=============================================================================== ! subroutine metablock_set_leaf(pmeta) implicit none type(block_meta), pointer, intent(inout) :: pmeta !------------------------------------------------------------------------------- ! ! return, if it is a leaf ! if (pmeta%leaf) return ! set the block's leaf flag ! pmeta%leaf = .true. ! increase the number of leafs ! nleafs = nleafs + 1 !------------------------------------------------------------------------------- ! end subroutine metablock_set_leaf ! !=============================================================================== ! ! subroutine METABLOCK_UNSET_LEAF: ! ------------------------------- ! ! Subroutine marks the meta block pointed by the input argument as non-leaf. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! !=============================================================================== ! subroutine metablock_unset_leaf(pmeta) implicit none type(block_meta), pointer, intent(inout) :: pmeta !------------------------------------------------------------------------------- ! ! return, if is not a leaf ! if (.not. pmeta%leaf) return ! unset the block's leaf flag ! pmeta%leaf = .false. ! decrease the number of leafs ! nleafs = nleafs - 1 !------------------------------------------------------------------------------- ! end subroutine metablock_unset_leaf ! !=============================================================================== ! ! subroutine METABLOCK_SET_UPDATE: ! ------------------------------- ! ! Subroutine marks the meta block pointed by the input argument to be updated. ! ! Arguments: ! ! pmeta - a pointer to the updated meta block; ! !=============================================================================== ! subroutine metablock_set_update(pmeta) implicit none type(block_meta), pointer, intent(inout) :: pmeta !------------------------------------------------------------------------------- ! pmeta%update = .true. !------------------------------------------------------------------------------- ! end subroutine metablock_set_update ! !=============================================================================== ! ! subroutine BUILD_LEAF_LIST: ! --------------------------- ! ! Subroutine builds the list of leaf blocks. ! ! Arguments: ! ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine build_leaf_list(status) use iso_fortran_env, only : error_unit implicit none integer, intent(out) :: status type(block_meta), pointer :: pmeta type(block_leaf), pointer :: pleaf character(len=*), parameter :: loc = 'BLOCKS::build_leaf_list()' !------------------------------------------------------------------------------- ! status = 0 ! if the list_leaf is associated, wipe it out ! if (associated(list_leaf)) then call wipe_leaf_list(status) if (status /= 0) then write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Cannot wipe the list of leafs!" status = 1 return end if end if ! associate pmeta with the first block on the meta block list ! pmeta => list_meta ! iterate over all meta blocks in the list ! do while (associated(pmeta)) ! check if the block is the leaf ! if (pmeta%leaf) then ! allocate new leaf structure ! allocate(pleaf, stat = status) ! check if allocation succeeded ! if (status == 0) then ! associate its pointers ! pleaf%next => list_leaf pleaf%meta => pmeta ! associate list_leaf to the allocated structure ! list_leaf => pleaf else write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Cannot allocate space for leaf!" status = 1 return end if end if ! leaf ! associate pmeta with the next meta block ! pmeta => pmeta%next end do ! meta blocks !------------------------------------------------------------------------------- ! end subroutine build_leaf_list ! !=============================================================================== ! ! subroutine BUILD_DATABLOCK_LIST: ! ------------------------------- ! ! Subroutine builds the list of data blocks data_blocks. ! ! Arguments: ! ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine build_datablock_list(status) use iso_fortran_env, only : error_unit implicit none integer, intent(out) :: status type(block_data), pointer :: pdata integer :: l character(len=*), parameter :: loc = 'BLOCKS::build_datablock_list()' !------------------------------------------------------------------------------- ! status = 0 if (allocated(data_blocks)) then deallocate(data_blocks, stat=status) if (status /= 0) & write(error_unit,"('[',a,']: ',a)") trim(loc), & "Could not deallocate the vector of data blocks!" end if allocate(data_blocks(dblocks), stat=status) if (status == 0) then l = 0 pdata => list_data do while (associated(pdata)) l = l + 1 data_blocks(l)%ptr => pdata pdata => pdata%next end do else write(error_unit,"('[',a,']: ',a)") trim(loc), & "Could not allocate the vector of data blocks!" end if !------------------------------------------------------------------------------- ! end subroutine build_datablock_list ! !=============================================================================== ! ! subroutine SET_NEIGHBORS_REFINE: ! ------------------------------- ! ! Subroutine marks all neighbors (including edge and corner ones) of ! the meta block pointed by the input argument to be refined if they ! fell under some certain conditions. ! ! Arguments: ! ! pmeta - a pointer to the refined meta block; ! !=============================================================================== ! subroutine set_neighbors_refine(pmeta) implicit none type(block_meta), pointer, intent(inout) :: pmeta procedure(reset_neighbors_update), pointer :: pprocedure !------------------------------------------------------------------------------- ! ! prepare the procedure pointer ! pprocedure => reset_neighbors_refinement ! iterate over all neighbors and coll pprocedure ! call iterate_over_neighbors(pmeta, pprocedure) !------------------------------------------------------------------------------- ! end subroutine set_neighbors_refine #ifdef DEBUG ! !=============================================================================== ! ! subroutine CHECK_NEIGHBORS: ! -------------------------- ! ! Subroutine iterates over all blocks and checks if the pointers to their ! neighbors are consistent. ! !=============================================================================== ! subroutine check_neighbors() implicit none type(block_meta), pointer :: pmeta !------------------------------------------------------------------------------- ! pmeta => list_meta do while(associated(pmeta)) ! check the block neighbors ! call check_block_neighbors(pmeta) pmeta => pmeta%next end do ! meta blocks !------------------------------------------------------------------------------- ! end subroutine check_neighbors #endif /* DEBUG */ ! !=============================================================================== !! !!*** 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; ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine insert_metablock_after(pprev, pmeta, status) use iso_fortran_env, only : error_unit implicit none type(block_meta), pointer, intent(in) :: pprev type(block_meta), pointer, intent(out) :: pmeta integer , intent(out) :: status character(len=*), parameter :: loc = 'BLOCKS::insert_metablock_after()' !------------------------------------------------------------------------------- ! ! allocate memory for the new meta block ! call allocate_metablock(pmeta, status) ! quit if block couldn't be allocated ! if (status /= 0) return ! 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 ! update the last_meta pointer if necessary ! if (associated(last_meta)) then if (pprev%id == last_meta%id) last_meta => pmeta else write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Argument pprev is associated but last_meta is not!" status = 1 return end if ! if pprev is null and list_meta is associated, there is something wrong ! else if (associated(list_meta)) then write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Argument pprev is null but list_meta is associated!" status = 1 return 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 ! increase the number of allocated meta blocks stored in the meta block list ! mblocks = mblocks + 1 !------------------------------------------------------------------------------- ! end subroutine insert_metablock_after ! !=============================================================================== ! ! subroutine WIPE_LEAF_LIST: ! -------------------------- ! ! Subroutine releases the memory used by the leaf block list. ! ! Arguments: ! ! status - the subroutine call status: 0 for success, otherwise failure; ! !=============================================================================== ! subroutine wipe_leaf_list(status) use iso_fortran_env, only : error_unit implicit none integer, intent(out) :: status type(block_leaf), pointer :: pleaf character(len=*), parameter :: loc = 'BLOCKS::wipe_leaf_list()' !------------------------------------------------------------------------------- ! status = 0 ! associate pleaf with the first pointer on the leaf list ! pleaf => list_leaf ! iterate over all leafs in the list ! do while (associated(pleaf)) ! associate the list pointer to the next leaf block ! list_leaf => pleaf%next ! deallocate the current leaf block ! deallocate(pleaf, stat = status) ! check if deallocation succeeded ! if (status /= 0) then write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Cannot deallocate leaf!" status = 1 return end if ! associate pleaf with the first pointer on the leaf list ! pleaf => list_leaf end do ! leaf list !------------------------------------------------------------------------------- ! end subroutine wipe_leaf_list ! !=============================================================================== ! ! function INCREASE_ID: ! -------------------- ! ! Function increases the last identification number by 1 and returns its ! value. ! !=============================================================================== ! function increase_id() result(id) implicit none integer(kind=4) :: id !------------------------------------------------------------------------------- ! last_id = last_id + 1 id = last_id !------------------------------------------------------------------------------- ! end function increase_id ! !=============================================================================== ! ! subroutine SET_NEIGHBORS_UPDATE: ! ------------------------------- ! ! Subroutine marks all neighbors (including edge and corner ones) of ! the meta block pointed by the input argument to be updated too. ! ! Arguments: ! ! pmeta - a pointer to the refined meta block; ! !=============================================================================== ! subroutine set_neighbors_update(pmeta) implicit none type(block_meta), pointer, intent(inout) :: pmeta procedure(reset_neighbors_update), pointer :: pprocedure !------------------------------------------------------------------------------- ! ! prepare the procedure pointer ! pprocedure => reset_neighbors_update ! iterate over all neighbors and coll pprocedure ! call iterate_over_neighbors(pmeta, pprocedure) !------------------------------------------------------------------------------- ! end subroutine set_neighbors_update ! !=============================================================================== ! ! subroutine RESET_NEIGHBORS_UPDATE: ! --------------------------------- ! ! Subroutine set the neighbor to be updated as well. ! ! Arguments: ! ! pmeta - a pointer to the refined meta block; ! pneigh - a pointer to the neighbor meta block; ! !=============================================================================== ! subroutine reset_neighbors_update(pmeta, pneigh) implicit none type(block_meta), pointer, intent(inout) :: pmeta, pneigh !------------------------------------------------------------------------------- ! if (associated(pmeta)) call metablock_set_update(pneigh) !------------------------------------------------------------------------------- ! end subroutine reset_neighbors_update ! !=============================================================================== ! ! subroutine RESET_NEIGHBORS_REFINEMENT: ! ------------------------------------- ! ! Subroutine checks the level of the neighbor block and depending on ! the refinement flags of both block resets it to the correct value. ! ! Arguments: ! ! pmeta - a pointer to the refined meta block; ! pneigh - a pointer to the neighbor meta block; ! !=============================================================================== ! subroutine reset_neighbors_refinement(pmeta, pneigh) implicit none type(block_meta), pointer, intent(inout) :: pmeta, pneigh !------------------------------------------------------------------------------- ! !=== conditions for blocks selected to be refined ! if (pmeta%refine == 1) then ! if the neighbor is set to be derefined, reset its flag (this applies to ! blocks at the current or lower level) ! pneigh%refine = max(0, pneigh%refine) ! if the neighbor is at lower level, always set it to be refined ! if (pneigh%level < pmeta%level) pneigh%refine = 1 end if ! refine = 1 !=== conditions for blocks which stay at the same level ! if (pmeta%refine == 0) then ! if the neighbor lays at lower level and is set to be derefined, cancel its ! derefinement ! if (pneigh%level < pmeta%level) pneigh%refine = max(0, pneigh%refine) end if ! refine = 0 !=== conditions for blocks which are selected to be derefined ! if (pmeta%refine == -1) then ! if the neighbor is at lower level and is set to be derefined, cancel its ! derefinement ! if (pneigh%level < pmeta%level) pneigh%refine = max(0, pneigh%refine) ! if a neighbor is set to be refined, cancel the derefinement of current block ! if (pneigh%refine == 1) pmeta%refine = 0 end if ! refine = -1 !------------------------------------------------------------------------------- ! end subroutine reset_neighbors_refinement ! !=============================================================================== ! ! subroutine ITERATE_OVER_NEIGHBORS: ! --------------------------------- ! ! Subroutine iterates over all neighbors of the meta block (including edge ! and corner ones) and executes a subroutine provided by the pointer. ! ! Arguments: ! ! pmeta - a pointer to the meta block which neighbors are iterated over; ! pproc - a pointer to the subroutine called with each pair (pmeta, pneigh); ! !=============================================================================== ! subroutine iterate_over_neighbors(pmeta, pprocedure) implicit none type(block_meta) , pointer, intent(inout) :: pmeta procedure(reset_neighbors_update), pointer, intent(in) :: pprocedure type(block_meta), pointer :: pneigh integer :: n, i, j #if NDIMS == 3 integer :: k #endif /* NDIMS == 3 */ !------------------------------------------------------------------------------- ! ! return if pmeta not associated ! if (.not. associated(pmeta)) return #if NDIMS == 2 ! iterate over edges and corners ! do j = 1, nsides do i = 1, nsides do n = 1, ndims ! associate pneigh with the face neighbor ! pneigh => pmeta%edges(i,j,n)%ptr ! call the procedure for the face neighbor ! if (associated(pneigh)) call pprocedure(pmeta, pneigh) end do ! n = 1, ndims ! associate pneigh with the face neighbor ! pneigh => pmeta%corners(i,j)%ptr ! call the procedure for the face neighbor ! if (associated(pneigh)) call pprocedure(pmeta, pneigh) end do ! i = 1, nsides end do ! j = 1, nsides #endif /* NDIMS == 2 */ #if NDIMS == 3 ! iterate over faces, edges, and corners ! do k = 1, nsides do j = 1, nsides do i = 1, nsides do n = 1, ndims ! associate pneigh with the face neighbor ! pneigh => pmeta%faces(i,j,k,n)%ptr ! call the procedure for the face neighbor ! if (associated(pneigh)) call pprocedure(pmeta, pneigh) ! associate pneigh with the face neighbor ! pneigh => pmeta%edges(i,j,k,n)%ptr ! call the procedure for the face neighbor ! if (associated(pneigh)) call pprocedure(pmeta, pneigh) end do ! n = 1, ndims ! associate pneigh with the face neighbor ! pneigh => pmeta%corners(i,j,k)%ptr ! call the procedure for the face neighbor ! if (associated(pneigh)) call pprocedure(pmeta, pneigh) end do ! i = 1, nsides end do ! j = 1, nsides end do ! k = 1, nsides #endif /* NDIMS == 3 */ !------------------------------------------------------------------------------- ! end subroutine iterate_over_neighbors #ifdef DEBUG ! !=============================================================================== ! ! subroutine CHECK_BLOCK_NEIGHBORS: ! -------------------------------- ! ! Subroutine iterates over all neighbor blocks and checks if their pointers ! are consistent, i.e. if their corresponding neighbor pointers point to ! the current block. ! ! Arguments: ! ! pmeta - a pointer to the meta block; ! !=============================================================================== ! subroutine check_block_neighbors(pmeta) use iso_fortran_env, only : error_unit implicit none type(block_meta), pointer, intent(in) :: pmeta type(block_meta), pointer :: pneigh, pself integer :: ip, ir, ic integer :: jp, jr, jc #if NDIMS == 3 integer :: kp, kr, kc #endif /* NDIMS == 2 */ character(len=*), parameter :: loc = 'BLOCKS::check_block_neighbors()' !------------------------------------------------------------------------------- ! ! return if the block is not leaf ! if (.not. pmeta%leaf) return #if NDIMS == 2 ! iterate over all corners ! do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip !--- check edges --- ! ! along X direction ! ! associate pneigh with the current corner ! pneigh => pmeta%edges(ip,jp,1)%ptr ! check if pneigh is associated ! if (associated(pneigh)) then ! if pneigh is on the same level ! if (pneigh%level == pmeta%level) then ! assiociate pself to the corresponding edge of the neighbor ! pself => pneigh%edges(ip,jr,1)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent same level neighbor edges!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, 1 write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Same level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, 1 end if ! pself associated end if ! pneigh and pmeta on the same level ! if pneigh is on the higher level level ! if (pneigh%level > pmeta%level) then ! iterate over all edges in the given direction and Y side ! do ic = 1, nsides ! assiociate pself to the corresponding edge of the neighbor ! pself => pneigh%edges(ic,jr,1)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent higher level neighbor edge!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, 1 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'self ', pself%id, ic, jr, 1 end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Higher level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ic, jr, 1 end if ! pself associated end do ! ic = 1, nsides end if ! pneigh on higher level end if ! pneigh associated ! along Y direction ! ! associate pneigh with the current corner ! pneigh => pmeta%edges(ip,jp,2)%ptr ! check if pneigh is associated ! if (associated(pneigh)) then ! if pneigh is on the same level ! if (pneigh%level == pmeta%level) then ! assiociate pself to the corresponding edge of the neighbor ! pself => pneigh%edges(ir,jp,2)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent same level neighbor edge!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, 2 write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Same level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, 2 end if ! pself associated end if ! pneigh and pmeta on the same level ! if pneigh is on the higher level level ! if (pneigh%level > pmeta%level) then ! iterate over all edges in the given direction and Y side ! do jc = 1, nsides ! assiociate pself to the corresponding edge of the neighbor ! pself => pneigh%edges(ir,jc,2)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent higher level neighbor edge!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'self ', pself%id, ir, jc, 2 end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Higher level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ir, jc, 2 end if ! pself associated end do ! jc = 1, nsides end if ! pneigh on higher level end if ! pneigh associated !--- check corners --- ! ! associate pneigh with the current corner ! pneigh => pmeta%corners(ip,jp)%ptr ! check if the neighbor is associated ! if (associated(pneigh)) then ! assiociate pself to the corresponding corner of the neighbor ! pself => pneigh%corners(ir,jr)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent neighbor corners!" write(error_unit,"(a6,' id: ',i8,' [ ',2(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp write(error_unit,"(a6,' id: ',i8,' [ ',2(i2,','),' ]')") & 'neigh', pneigh%id, ir, jr write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Neighbor's corner not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',2(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp write(error_unit,"(a6,' id: ',i8,' [ ',2(i2,','),' ]')") & 'neigh', pneigh%id, ir, jr end if ! pself associated end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides #endif /* NDIMS == 2 */ #if NDIMS == 3 ! iterate over all corners ! do kp = 1, nsides kr = 3 - kp do jp = 1, nsides jr = 3 - jp do ip = 1, nsides ir = 3 - ip !--- check faces --- ! ! along X direction ! ! associate pneigh with the current face ! pneigh => pmeta%faces(ip,jp,kp,1)%ptr ! check if pneigh is associated ! if (associated(pneigh)) then ! if pneigh is on the same level ! if (pneigh%level == pmeta%level) then ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%faces(ir,jp,kp,1)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent same level neighbor faces!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, kp, 1 write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Same level neighbor's face not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, kp, 1 end if ! pself associated end if ! pneigh and pmeta on the same level ! if pneigh is on higher level ! if (pneigh%level > pmeta%level) then ! iterate over all neighbor faces ! do kc = 1, nsides do jc = 1, nsides ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%faces(ir,jc,kc,1)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent higher level neighbor face!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'self ', pself%id, ir, jc, kc, 1 end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Higher level neighbor's face not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jc, kc, 1 end if ! pself associated end do ! jc = 1, nsides end do ! kc = 1, nsides end if ! pneigh on higher level end if ! pneigh associated ! along Y direction ! ! associate pneigh with the current face ! pneigh => pmeta%faces(ip,jp,kp,2)%ptr ! check if pneigh is associated ! if (associated(pneigh)) then ! if pneigh is on the same level ! if (pneigh%level == pmeta%level) then ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%faces(ip,jr,kp,2)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent same level neighbor faces!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, kp, 2 write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Same level neighbor's face not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, kp, 2 end if ! pself associated end if ! pneigh and pmeta on the same level ! if pneigh is on higher level ! if (pneigh%level > pmeta%level) then ! iterate over all neighbor faces ! do kc = 1, nsides do ic = 1, nsides ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%faces(ic,jr,kc,2)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent higher level neighbor face!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'self ', pself%id, ic, jr, kc, 2 end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Higher level neighbor's face not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ic, jr, kc, 2 end if ! pself associated end do ! ic = 1, nsides end do ! kc = 1, nsides end if ! pneigh on higher level end if ! pneigh associated ! along Z direction ! ! associate pneigh with the current face ! pneigh => pmeta%faces(ip,jp,kp,3)%ptr ! check if pneigh is associated ! if (associated(pneigh)) then ! if pneigh is on the same level ! if (pneigh%level == pmeta%level) then ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%faces(ip,jp,kr,3)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent same level neighbor faces!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jp, kr, 3 write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Same level neighbor's face not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jp, kr, 3 end if ! pself associated end if ! pneigh and pmeta on the same level ! if pneigh is on higher level ! if (pneigh%level > pmeta%level) then ! iterate over all neighbor faces ! do jc = 1, nsides do ic = 1, nsides ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%faces(ic,jc,kr,3)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent higher level neighbor face!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jp, kr, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'self ', pself%id, ic, jc, kr, 3 end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Higher level neighbor's face not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ic, jc, kr, 3 end if ! pself associated end do ! ic = 1, nsides end do ! jc = 1, nsides end if ! pneigh on higher level end if ! pneigh associated !--- check edges --- ! ! along X direction ! ! associate pneigh with the current edge ! pneigh => pmeta%edges(ip,jp,kp,1)%ptr ! check if pneigh is associated ! if (associated(pneigh)) then ! if pneigh is on the same level ! if (pneigh%level == pmeta%level) then ! assiociate pself to the corresponding edge of the neighbor ! pself => pneigh%edges(ip,jr,kr,1)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent same level neighbor edges!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, kr, 1 write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Same level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, kr, 1 end if ! pself associated end if ! pneigh and pmeta on the same level ! if pneigh is on higher level ! if (pneigh%level > pmeta%level) then ! iterate over all neighbor edges ! do ic = 1, nsides ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%edges(ic,jr,kr,1)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent higher level neighbor edge!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ip, jr, kr, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'self ', pself%id, ic, jr, kr, 1 end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Higher level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 1 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ic, jr, kr, 1 end if ! pself associated end do ! ic = 1, nsides end if ! pneigh on higher level end if ! pneigh associated ! along Y direction ! ! associate pneigh with the current edge ! pneigh => pmeta%edges(ip,jp,kp,2)%ptr ! check if pneigh is associated ! if (associated(pneigh)) then ! if pneigh is on the same level ! if (pneigh%level == pmeta%level) then ! assiociate pself to the corresponding edge of the neighbor ! pself => pneigh%edges(ir,jp,kr,2)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent same level neighbor edges!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, kr, 2 write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Same level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, kr, 2 end if ! pself associated end if ! pneigh and pmeta on the same level ! if pneigh is on higher level ! if (pneigh%level > pmeta%level) then ! iterate over all neighbor edges ! do jc = 1, nsides ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%edges(ir,jc,kr,2)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent higher level neighbor edge!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jp, kr, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'self ', pself%id, ir, jc, kr, 2 end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Higher level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 2 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jc, kr, 2 end if ! pself associated end do ! jc = 1, nsides end if ! pneigh on higher level end if ! pneigh associated ! along Z direction ! ! associate pneigh with the current edge ! pneigh => pmeta%edges(ip,jp,kp,3)%ptr ! check if pneigh is associated ! if (associated(pneigh)) then ! if pneigh is on the same level ! if (pneigh%level == pmeta%level) then ! assiociate pself to the corresponding edge of the neighbor ! pself => pneigh%edges(ir,jr,kp,3)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent same level neighbor edges!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jr, kp, 3 write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Same level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jr, kp, 3 end if ! pself associated end if ! pneigh and pmeta on the same level ! if pneigh is on higher level ! if (pneigh%level > pmeta%level) then ! iterate over all neighbor edges ! do kc = 1, nsides ! assiociate pself to the corresponding face of the neighbor ! pself => pneigh%edges(ir,jr,kc,3)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent higher level neighbor edge!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jr, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'self ', pself%id, ir, jr, kc, 3 end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Higher level neighbor's edge not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'meta ', pmeta%id , ip, jp, kp, 3 write(error_unit,"(a6,' id: ',i8,' [ ',4(i2,','),' ]')") & 'neigh', pneigh%id, ir, jr, kc, 3 end if ! pself associated end do ! kc = 1, nsides end if ! pneigh on higher level end if ! pneigh associated !--- check corners --- ! ! associate pneigh with the current corner ! pneigh => pmeta%corners(ip,jp,kp)%ptr ! check if the neighbor is associated ! if (associated(pneigh)) then ! assiociate pself to the corresponding corner of the neighbor ! pself => pneigh%corners(ir,jr,kr)%ptr ! check if pself is associated ! if (associated(pself)) then ! check if pself is the same as pmeta ! if (pmeta%id /= pself%id) then ! print warning, since the blocks differ ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Inconsistent neighbor corners!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp, kp write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ir, jr, kr write(error_unit,"(a6,' id: ',i8)") 'self ', pself%id end if ! %id fields don't match else ! pself associated ! print warning, since the pointer should be associated ! write(error_unit,"('[',a,']: ',a)") trim(loc) & , "Neighbor's corner not associated!" write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'meta ', pmeta%id, ip, jp, kp write(error_unit,"(a6,' id: ',i8,' [ ',3(i2,','),' ]')") & 'neigh', pneigh%id, ir, jr, kr end if ! pself associated end if ! pneigh associated end do ! ip = 1, nsides end do ! jp = 1, nsides end do ! kp = 1, nsides #endif /* NDIMS == 3 */ !------------------------------------------------------------------------------- ! end subroutine check_block_neighbors #endif /* DEBUG */ !=============================================================================== ! end module