2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
2011-05-05 18:37:53 -03:00
|
|
|
!! module: MESH - handling adaptive mesh structure
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
2011-05-05 18:37:53 -03:00
|
|
|
!! Copyright (C) 2008-2011 Grzegorz Kowal <grzegorz@amuncode.org>
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
2011-05-05 18:37:53 -03:00
|
|
|
!! This file is part of the AMUN code.
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
2011-04-29 11:21:30 -03:00
|
|
|
!! This program is free software; you can redistribute it and/or
|
|
|
|
!! modify it under the terms of the GNU General Public License
|
|
|
|
!! as published by the Free Software Foundation; either version 2
|
|
|
|
!! of the License, or (at your option) any later version.
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
2011-04-29 11:21:30 -03:00
|
|
|
!! This program is distributed in the hope that it will be useful,
|
2008-11-11 16:12:26 -06:00
|
|
|
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
!! GNU General Public License for more details.
|
|
|
|
!!
|
|
|
|
!! You should have received a copy of the GNU General Public License
|
2011-04-29 11:21:30 -03:00
|
|
|
!! along with this program; if not, write to the Free Software Foundation,
|
|
|
|
!! Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-11-11 16:12:26 -06:00
|
|
|
!!
|
|
|
|
!
|
|
|
|
module mesh
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2008-12-09 14:51:33 -06:00
|
|
|
! minimum grid step
|
|
|
|
!
|
|
|
|
real, save :: dx_min = 1.0
|
|
|
|
|
2009-05-18 22:46:19 +02:00
|
|
|
! spatial coordinates for all levels of refinements
|
2008-12-08 20:03:01 -06:00
|
|
|
!
|
|
|
|
real, dimension(:,:), allocatable, save :: ax , ay , az
|
|
|
|
real, dimension(: ), allocatable, save :: adx , ady , adz
|
|
|
|
real, dimension(: ), allocatable, save :: adxi, adyi, adzi
|
2011-03-02 15:01:08 -03:00
|
|
|
real, dimension(: ), allocatable, save :: advol
|
2008-12-08 20:03:01 -06:00
|
|
|
|
2011-04-28 12:04:22 -03:00
|
|
|
! maximum number of covering blocks
|
|
|
|
!
|
|
|
|
integer , save :: tblocks = 1
|
|
|
|
integer,dimension(NDIMS), save :: effres = 1
|
|
|
|
|
2011-05-05 16:51:35 -03:00
|
|
|
! the effective resolution for all levels
|
|
|
|
!
|
|
|
|
integer(kind=4), dimension(:,:), allocatable, save :: res
|
|
|
|
|
2011-04-28 12:04:22 -03:00
|
|
|
! log file for the mesh statistics
|
|
|
|
!
|
|
|
|
character(len=32), save :: fname = "mesh.log"
|
|
|
|
integer , save :: funit = 11
|
|
|
|
|
2008-11-11 16:12:26 -06:00
|
|
|
contains
|
|
|
|
!
|
2008-12-13 15:08:18 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2011-04-27 22:39:25 -03:00
|
|
|
! init_mesh: subroutine initializes mesh module and coordinate variables
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-04-28 12:04:22 -03:00
|
|
|
subroutine init_mesh(flag)
|
2011-04-27 22:39:25 -03:00
|
|
|
|
2011-05-05 18:25:24 -03:00
|
|
|
use blocks , only : datablock_set_dims
|
2011-05-05 16:51:35 -03:00
|
|
|
use config , only : maxlev, in, jn, kn, im, jm, km, ncells, rdims, ng &
|
2011-05-05 15:56:38 -03:00
|
|
|
, xmin, xmax, ymin, ymax, zmin, zmax
|
|
|
|
use mpitools , only : is_master, ncpus
|
|
|
|
use timer , only : start_timer, stop_timer
|
|
|
|
use variables, only : nqt, nvr
|
2011-04-27 22:39:25 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2011-04-28 12:04:22 -03:00
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
logical, intent(in) :: flag
|
|
|
|
|
2011-04-27 22:39:25 -03:00
|
|
|
! local variables
|
|
|
|
!
|
2011-04-28 12:04:22 -03:00
|
|
|
integer(kind=4) :: i, j, k, l, n
|
|
|
|
logical :: info
|
2011-04-27 22:39:25 -03:00
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
integer(kind=4), dimension(3) :: bm, rm, dm
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-05-02 23:35:48 -03:00
|
|
|
! start the mesh timer
|
|
|
|
!
|
|
|
|
call start_timer(5)
|
|
|
|
|
2011-05-05 15:56:38 -03:00
|
|
|
! set data block dimensions
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
call datablock_set_dims(nqt, nvr, im, jm, km)
|
2011-05-05 15:56:38 -03:00
|
|
|
|
2011-04-27 22:39:25 -03:00
|
|
|
! allocating space for coordinate variables
|
|
|
|
!
|
|
|
|
allocate(ax (maxlev, im))
|
|
|
|
allocate(ay (maxlev, jm))
|
|
|
|
allocate(az (maxlev, km))
|
|
|
|
allocate(adx (maxlev))
|
|
|
|
allocate(ady (maxlev))
|
|
|
|
allocate(adz (maxlev))
|
|
|
|
allocate(adxi (maxlev))
|
|
|
|
allocate(adyi (maxlev))
|
|
|
|
allocate(adzi (maxlev))
|
|
|
|
allocate(advol(maxlev))
|
2011-05-05 16:51:35 -03:00
|
|
|
allocate(res (maxlev, 3))
|
2011-04-27 22:39:25 -03:00
|
|
|
|
|
|
|
! reset the coordinate variables
|
|
|
|
!
|
|
|
|
ax(:,:) = 0.0d0
|
|
|
|
ay(:,:) = 0.0d0
|
|
|
|
az(:,:) = 0.0d0
|
|
|
|
adx(:) = 1.0d0
|
|
|
|
ady(:) = 1.0d0
|
|
|
|
adz(:) = 1.0d0
|
|
|
|
adxi(:) = 1.0d0
|
|
|
|
adyi(:) = 1.0d0
|
|
|
|
adzi(:) = 1.0d0
|
|
|
|
advol(:) = 1.0d0
|
2011-05-05 16:51:35 -03:00
|
|
|
res(:,:) = 1
|
2011-04-27 22:39:25 -03:00
|
|
|
|
|
|
|
! generating coordinates for all levels
|
|
|
|
!
|
|
|
|
do l = 1, maxlev
|
2011-05-05 16:51:35 -03:00
|
|
|
|
2011-04-27 22:39:25 -03:00
|
|
|
n = ncells * 2**(l - 1)
|
|
|
|
|
2011-05-05 16:51:35 -03:00
|
|
|
! calculate the effective resolution at each level
|
|
|
|
!
|
|
|
|
res(l,1) = in * 2**(maxlev - l)
|
|
|
|
res(l,2) = jn * 2**(maxlev - l)
|
|
|
|
#if NDIMS == 3
|
|
|
|
res(l,3) = kn * 2**(maxlev - l)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
2011-04-27 22:39:25 -03:00
|
|
|
adx (l) = (xmax - xmin) / (rdims(1) * n)
|
|
|
|
ady (l) = (ymax - ymin) / (rdims(2) * n)
|
|
|
|
#if NDIMS == 3
|
|
|
|
adz (l) = (zmax - zmin) / (rdims(3) * n)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
ax(l,:) = ((/(i, i = 1, im)/) - ng - 0.5d0) * adx(l)
|
|
|
|
ay(l,:) = ((/(j, j = 1, jm)/) - ng - 0.5d0) * ady(l)
|
|
|
|
#if NDIMS == 3
|
|
|
|
az(l,:) = ((/(k, k = 1, km)/) - ng - 0.5d0) * adz(l)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
adxi(l) = 1.0d0 / adx(l)
|
|
|
|
adyi(l) = 1.0d0 / ady(l)
|
|
|
|
#if NDIMS == 3
|
|
|
|
adzi(l) = 1.0d0 / adz(l)
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
advol(l) = adx(l) * ady(l) * adz(l)
|
|
|
|
end do
|
|
|
|
|
|
|
|
! get the minimum grid step
|
|
|
|
!
|
|
|
|
dx_min = 0.5 * min(adx(maxlev), ady(maxlev), adz(maxlev))
|
|
|
|
|
|
|
|
! print general information about resolutions
|
|
|
|
!
|
|
|
|
if (is_master()) then
|
|
|
|
|
|
|
|
bm( : ) = 1
|
|
|
|
rm( : ) = 1
|
|
|
|
dm( : ) = 1
|
|
|
|
|
|
|
|
bm(1:NDIMS) = rdims(1:NDIMS) * ncells
|
2011-05-05 16:51:35 -03:00
|
|
|
rm(1:NDIMS) = rdims(1:NDIMS) * res(1,1:NDIMS)
|
2011-04-27 22:39:25 -03:00
|
|
|
dm(1:NDIMS) = rm(1:NDIMS) / ncells
|
|
|
|
|
2011-04-28 12:04:22 -03:00
|
|
|
effres(1:NDIMS) = rm(1:NDIMS)
|
|
|
|
tblocks = product(dm(:))
|
|
|
|
|
2011-04-27 22:39:25 -03:00
|
|
|
write(*,*)
|
|
|
|
write(*,"(1x,a)" ) "Geometry:"
|
|
|
|
write(*,"(4x,a, 1x,i6)" ) "refinement to level =", maxlev
|
|
|
|
write(*,"(4x,a,3(1x,i6))") "base configuration =", rdims(1:NDIMS)
|
|
|
|
write(*,"(4x,a,3(1x,i6))") "top level blocks =", dm(1:NDIMS)
|
2011-04-28 12:04:22 -03:00
|
|
|
write(*,"(4x,a, 1x,i6)" ) "maxium cover blocks =", tblocks
|
2011-04-27 22:39:25 -03:00
|
|
|
write(*,"(4x,a,3(1x,i6))") "base resolution =", bm(1:NDIMS)
|
|
|
|
write(*,"(4x,a,3(1x,i6))") "effective resolution =", rm(1:NDIMS)
|
|
|
|
|
2011-04-28 12:04:22 -03:00
|
|
|
! prepare the file for logging mesh statistics; only the master process handles
|
|
|
|
! this part
|
|
|
|
!
|
|
|
|
! check if the mesh statistics file exists
|
|
|
|
!
|
|
|
|
inquire(file = fname, exist = info)
|
|
|
|
|
|
|
|
! if flag is set to .true. or the original mesh statistics file does not exist,
|
|
|
|
! create a new one and write the header in it, otherwise open the original file
|
|
|
|
! and move the writing position to the end to allow for appending
|
|
|
|
!
|
|
|
|
if (flag .or. .not. info) then
|
|
|
|
|
|
|
|
! create a new mesh statistics file
|
|
|
|
!
|
|
|
|
open(unit=funit, file=fname, form='formatted', status='replace')
|
|
|
|
|
|
|
|
! write the mesh statistics file header
|
|
|
|
!
|
|
|
|
write(funit,"('#')")
|
|
|
|
write(funit,"('#',4x,'step',5x,'time',11x,'leaf',4x,'meta'" // &
|
|
|
|
"6x,'coverage',8x,'AMR',11x,'block distribution')")
|
|
|
|
write(funit,"('#',27x,'blocks',2x,'blocks',4x,'efficiency'" // &
|
|
|
|
",6x,'efficiency',$)")
|
|
|
|
write(funit,"(1x,a12,$)") 'level = 1'
|
|
|
|
do l = 2, maxlev
|
|
|
|
write(funit,"(2x,i6,$)") l
|
|
|
|
end do
|
|
|
|
#ifdef MPI
|
|
|
|
write(funit,"(1x,a10,$)") 'cpu = 1'
|
|
|
|
do n = 2, ncpus
|
|
|
|
write(funit,"(2x,i6,$)") n
|
|
|
|
end do
|
|
|
|
#endif /* MPI */
|
|
|
|
write(funit,"('' )")
|
|
|
|
write(funit,"('#')")
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! open the mesh statistics file and set the position at the end of file
|
|
|
|
!
|
|
|
|
open(unit=funit, file=fname, form='formatted', position='append')
|
|
|
|
|
|
|
|
! write a marker that the job has been restarted from here
|
|
|
|
!
|
|
|
|
write(funit,"('#',1x,a)") "job restarted from this point"
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
end if ! master
|
2011-05-02 23:35:48 -03:00
|
|
|
|
|
|
|
! stop the mesh timer
|
2011-04-27 22:39:25 -03:00
|
|
|
!
|
2011-05-02 23:35:48 -03:00
|
|
|
call stop_timer(5)
|
|
|
|
|
2011-04-27 22:39:25 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine init_mesh
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-04-27 22:35:38 -03:00
|
|
|
! generate_mesh: subroutine generate the initial block structure by creating
|
|
|
|
! blocks according to the geometry, initial problem and
|
|
|
|
! refinement criterium
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2008-12-13 15:08:18 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2011-04-27 22:35:38 -03:00
|
|
|
subroutine generate_mesh()
|
2008-11-11 16:12:26 -06:00
|
|
|
|
2011-04-15 00:46:19 +02:00
|
|
|
use config , only : maxlev, rdims
|
2011-04-14 00:19:45 +02:00
|
|
|
use blocks , only : block_meta, block_data, list_meta, list_data
|
|
|
|
use blocks , only : refine_block, deallocate_datablock
|
2011-05-05 18:25:24 -03:00
|
|
|
use blocks , only : nchild, nsides, nfaces
|
|
|
|
use blocks , only : get_mblocks, get_nleafs
|
2008-12-28 13:09:14 -06:00
|
|
|
use error , only : print_info, print_error
|
|
|
|
use mpitools, only : is_master, ncpu, ncpus
|
2008-12-22 16:08:51 -06:00
|
|
|
use problem , only : init_domain, init_problem, check_ref
|
2011-05-02 23:35:48 -03:00
|
|
|
use timer , only : start_timer, stop_timer
|
2008-11-11 16:12:26 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2009-09-14 18:19:27 -03:00
|
|
|
! local pointers
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2009-09-14 19:15:21 -03:00
|
|
|
type(block_meta), pointer :: pmeta_block, pneigh, pnext
|
2009-09-11 21:52:18 -03:00
|
|
|
type(block_data), pointer :: pdata_block
|
2009-09-14 18:19:27 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2011-04-28 12:18:55 -03:00
|
|
|
integer(kind=4) :: i, j, k, l, n
|
2010-12-14 11:19:11 -02:00
|
|
|
#ifdef MPI
|
|
|
|
integer(kind=4), dimension(0:ncpus-1) :: lb
|
|
|
|
#endif /* MPI */
|
2008-11-11 16:12:26 -06:00
|
|
|
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2011-05-02 23:35:48 -03:00
|
|
|
! start the mesh timer
|
|
|
|
!
|
|
|
|
call start_timer(5)
|
|
|
|
|
2008-12-22 15:34:02 -06:00
|
|
|
! allocate the initial structure of blocks according to the problem
|
|
|
|
!
|
2011-05-05 16:51:35 -03:00
|
|
|
call init_domain(res(1,:))
|
2008-11-11 16:12:26 -06:00
|
|
|
|
|
|
|
! at this point we assume, that the initial structure of blocks
|
|
|
|
! according to the defined geometry is already created; no refinement
|
2009-05-18 22:46:19 +02:00
|
|
|
! is done yet; we fill out the coarse blocks with the initial condition
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2011-04-25 13:31:41 -03:00
|
|
|
if (is_master()) then
|
|
|
|
write(*,*)
|
|
|
|
write(*,"(1x,a)" ) "Generating the initial mesh:"
|
2011-02-16 12:08:50 -02:00
|
|
|
write(*,"(4x,a,$)") "generating level = "
|
2011-04-25 13:31:41 -03:00
|
|
|
end if
|
2008-12-07 12:56:00 -06:00
|
|
|
|
2008-12-28 13:09:14 -06:00
|
|
|
l = 1
|
2009-05-18 22:46:19 +02:00
|
|
|
do while (l .le. maxlev)
|
2008-12-28 13:09:14 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! print the level currently processed
|
2008-12-07 12:56:00 -06:00
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
if (is_master()) &
|
|
|
|
write(*,"(1x,i2,$)") l
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! iterate over all data blocks at the current level and initialize the problem
|
2008-11-29 22:19:02 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pdata_block => list_data
|
|
|
|
do while (associated(pdata_block))
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! set the initial conditions at the current block
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
if (pdata_block%meta%level .le. l) &
|
|
|
|
call init_problem(pdata_block)
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2008-12-06 20:11:36 -06:00
|
|
|
! assign pointer to the next block
|
2008-12-05 14:04:12 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pdata_block => pdata_block%next
|
2008-12-06 20:11:36 -06:00
|
|
|
end do
|
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! at the maximum level we only initialize the problem (without checking the
|
|
|
|
! refinement criterion)
|
2009-05-18 22:46:19 +02:00
|
|
|
!
|
|
|
|
if (l .lt. maxlev) then
|
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! iterate over all data blocks at the current level and check the refinement
|
2009-05-18 22:46:19 +02:00
|
|
|
! criterion; do not allow for derefinement
|
2008-12-06 20:11:36 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pdata_block => list_data
|
|
|
|
do while (associated(pdata_block))
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-09-29 16:17:10 -03:00
|
|
|
if (pdata_block%meta%level .eq. l) then
|
2009-09-11 21:52:18 -03:00
|
|
|
pdata_block%meta%refine = max(0, check_ref(pdata_block))
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-29 16:17:10 -03:00
|
|
|
! if there is only one block, and it is set not to be refined, refine it anyway
|
|
|
|
! because the resolution for the problem initiation may be too small
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
if (get_mblocks() .eq. 1 .and. l .eq. 1) &
|
2009-09-29 16:17:10 -03:00
|
|
|
pdata_block%meta%refine = 1
|
|
|
|
end if
|
|
|
|
|
2009-05-18 22:46:19 +02:00
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pdata_block => pdata_block%next
|
2009-05-18 22:46:19 +02:00
|
|
|
end do
|
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! walk through all levels down from the current level and check if select all
|
|
|
|
! neighbors for the refinement if they are at lower level; there is no need for
|
|
|
|
! checking the blocks at the lowest level;
|
2009-05-18 22:46:19 +02:00
|
|
|
!
|
|
|
|
do n = l, 2, -1
|
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! iterate over all meta blocks at the level n and if the current block is
|
2009-05-18 22:46:19 +02:00
|
|
|
! selected for the refinement and its neighbors are at lower levels select them
|
|
|
|
! for refinement too;
|
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pmeta_block => list_meta
|
|
|
|
do while (associated(pmeta_block))
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! check if the current block is at the level n, is a leaf, and selected for
|
|
|
|
! refinement
|
2008-12-06 20:11:36 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
if (pmeta_block%level .eq. n) then
|
|
|
|
if (pmeta_block%leaf) then
|
|
|
|
if (pmeta_block%refine .eq. 1) then
|
2008-12-06 20:11:36 -06:00
|
|
|
|
|
|
|
! iterate over all neighbors
|
|
|
|
!
|
2011-04-15 00:46:19 +02:00
|
|
|
do i = 1, NDIMS
|
2009-09-11 21:52:18 -03:00
|
|
|
do j = 1, nsides
|
|
|
|
do k = 1, nfaces
|
2008-12-28 13:09:14 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! assign pointer to the neighbor
|
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pneigh => pmeta_block%neigh(i,j,k)%ptr
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2008-12-28 13:09:14 -06:00
|
|
|
! check if the neighbor is associated
|
|
|
|
!
|
2009-05-18 22:46:19 +02:00
|
|
|
if (associated(pneigh)) then
|
2009-09-22 19:32:57 -03:00
|
|
|
|
|
|
|
! check if the neighbor is a leaf, if not something wrong is going on
|
|
|
|
!
|
2009-05-18 22:46:19 +02:00
|
|
|
if (pneigh%leaf) then
|
|
|
|
|
2009-09-11 21:52:18 -03:00
|
|
|
! if the neighbor has lower level, select it to be refined too
|
2008-12-05 14:43:43 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
if (pneigh%level .lt. pmeta_block%level) &
|
2009-05-18 22:46:19 +02:00
|
|
|
pneigh%refine = 1
|
2009-09-22 19:32:57 -03:00
|
|
|
|
2009-05-18 22:46:19 +02:00
|
|
|
else
|
|
|
|
call print_error("mesh::init_mesh", "Neighbor is not a leaf!")
|
2009-09-22 19:32:57 -03:00
|
|
|
end if
|
|
|
|
end if
|
2008-12-05 14:43:43 -06:00
|
|
|
|
2009-05-18 22:46:19 +02:00
|
|
|
end do
|
2008-12-28 13:09:14 -06:00
|
|
|
end do
|
|
|
|
end do
|
2008-12-05 14:43:43 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
end if
|
|
|
|
end if
|
|
|
|
end if
|
2008-12-06 20:11:36 -06:00
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pmeta_block => pmeta_block%next
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2009-05-18 22:46:19 +02:00
|
|
|
end do
|
2008-12-06 20:11:36 -06:00
|
|
|
end do
|
2008-12-05 14:43:43 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
!! refine all selected blocks starting from the lowest level
|
|
|
|
!!
|
|
|
|
! walk through the levels starting from the lowest to the current level
|
2008-12-05 14:43:43 -06:00
|
|
|
!
|
2009-05-18 22:46:19 +02:00
|
|
|
do n = 1, l
|
2008-12-05 14:43:43 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! iterate over all meta blocks
|
2008-12-05 14:43:43 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pmeta_block => list_meta
|
|
|
|
do while (associated(pmeta_block))
|
2008-12-05 14:43:43 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! check if the current block is at the level n and, if it is selected for
|
|
|
|
! refinement and if so, perform the refinement on this block
|
|
|
|
!
|
|
|
|
if (pmeta_block%level .eq. n .and. pmeta_block%refine .eq. 1) then
|
|
|
|
|
|
|
|
! perform the refinement
|
2008-12-06 20:11:36 -06:00
|
|
|
!
|
2011-05-05 16:51:35 -03:00
|
|
|
call refine_block(pmeta_block, res(pmeta_block%level + 1,:), .true.)
|
2008-12-05 14:43:43 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
end if
|
|
|
|
|
2008-12-28 13:09:14 -06:00
|
|
|
! assign pointer to the next block
|
2008-12-06 20:11:36 -06:00
|
|
|
!
|
2009-09-11 21:52:18 -03:00
|
|
|
pmeta_block => pmeta_block%next
|
2009-05-18 22:46:19 +02:00
|
|
|
end do
|
2008-12-28 13:09:14 -06:00
|
|
|
end do
|
2009-09-22 19:32:57 -03:00
|
|
|
end if
|
2008-12-06 20:11:36 -06:00
|
|
|
|
2008-12-28 13:09:14 -06:00
|
|
|
l = l + 1
|
|
|
|
end do
|
|
|
|
|
2009-09-26 00:39:41 -03:00
|
|
|
! deallocate data blocks of non leafs
|
|
|
|
!
|
|
|
|
pmeta_block => list_meta
|
|
|
|
do while (associated(pmeta_block))
|
|
|
|
|
|
|
|
if (.not. pmeta_block%leaf) &
|
|
|
|
call deallocate_datablock(pmeta_block%data)
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pmeta_block => pmeta_block%next
|
|
|
|
end do
|
|
|
|
|
2008-12-28 13:09:14 -06:00
|
|
|
#ifdef MPI
|
2009-09-25 20:35:38 -03:00
|
|
|
! divide blocks between all processes, use the number of data blocks to do this,
|
|
|
|
! but keep blocks from the top level which have the same parent packed together
|
2009-09-14 19:15:21 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
l = mod(get_nleafs(), ncpus) - 1
|
|
|
|
lb( : ) = get_nleafs() / ncpus
|
2010-12-14 11:19:11 -02:00
|
|
|
lb(0:l) = lb(0:l) + 1
|
|
|
|
|
|
|
|
! reset the processor and block numbers
|
2009-09-22 19:32:57 -03:00
|
|
|
n = 0
|
2010-12-14 11:19:11 -02:00
|
|
|
l = 0
|
2009-09-25 20:35:38 -03:00
|
|
|
|
2009-09-14 19:15:21 -03:00
|
|
|
pmeta_block => list_meta
|
|
|
|
do while (associated(pmeta_block))
|
|
|
|
|
|
|
|
! assign the cpu to the current block
|
|
|
|
!
|
2009-09-25 20:35:38 -03:00
|
|
|
pmeta_block%cpu = n
|
2009-09-14 19:15:21 -03:00
|
|
|
|
2009-09-25 20:35:38 -03:00
|
|
|
! increase the number of blocks on the current process; if it exceeds the
|
|
|
|
! allowed number reset the counter and increase the processor number
|
2009-09-14 19:15:21 -03:00
|
|
|
!
|
2009-09-26 00:39:41 -03:00
|
|
|
if (pmeta_block%leaf) then
|
2010-12-14 11:19:11 -02:00
|
|
|
l = l + 1
|
|
|
|
if (l .ge. lb(n)) then
|
2009-09-26 00:39:41 -03:00
|
|
|
n = min(ncpus - 1, n + 1)
|
2010-12-14 11:19:11 -02:00
|
|
|
l = 0
|
2009-09-25 20:35:38 -03:00
|
|
|
end if
|
|
|
|
end if
|
2009-09-14 19:15:21 -03:00
|
|
|
|
2009-09-25 20:35:38 -03:00
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pmeta_block => pmeta_block%next
|
2009-09-14 19:15:21 -03:00
|
|
|
end do
|
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! remove all data blocks which do not belong to the current process
|
2009-09-14 19:15:21 -03:00
|
|
|
!
|
|
|
|
pmeta_block => list_meta
|
|
|
|
do while (associated(pmeta_block))
|
|
|
|
pnext => pmeta_block%next
|
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! if the current block belongs to another process and its data field is
|
|
|
|
! associated, deallocate its data field
|
|
|
|
!
|
|
|
|
if (pmeta_block%cpu .ne. ncpu .and. associated(pmeta_block%data)) &
|
2009-09-14 19:15:21 -03:00
|
|
|
call deallocate_datablock(pmeta_block%data)
|
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
2009-09-14 19:15:21 -03:00
|
|
|
pmeta_block => pnext
|
|
|
|
end do
|
2008-12-28 13:09:14 -06:00
|
|
|
#endif /* MPI */
|
2008-11-29 22:19:02 -06:00
|
|
|
|
2011-04-28 12:18:55 -03:00
|
|
|
! go to a new line after generating levels
|
2008-12-07 12:56:00 -06:00
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
if (is_master()) then
|
2009-09-22 19:32:57 -03:00
|
|
|
write(*,*)
|
|
|
|
end if
|
2008-11-11 16:12:26 -06:00
|
|
|
|
2011-05-02 23:35:48 -03:00
|
|
|
! stop the mesh timer
|
|
|
|
!
|
|
|
|
call stop_timer(5)
|
|
|
|
|
2011-04-15 00:46:19 +02:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2011-04-27 22:35:38 -03:00
|
|
|
end subroutine generate_mesh
|
2011-04-15 00:46:19 +02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2008-12-13 15:08:18 -06:00
|
|
|
! update_mesh: subroutine check the refinement criterion for each block,
|
|
|
|
! refines or derefines it if necessary, and restricts or
|
|
|
|
! prolongates all data to the newly created blocks
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-02-27 13:49:26 -03:00
|
|
|
subroutine update_mesh()
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2010-12-01 09:25:30 -02:00
|
|
|
use config , only : maxlev, im, jm, km
|
|
|
|
use blocks , only : block_meta, block_data, list_meta, list_data &
|
2011-05-05 18:25:24 -03:00
|
|
|
, nchild, ndims, nsides, nfaces &
|
2010-12-01 09:25:30 -02:00
|
|
|
, refine_block, derefine_block, append_datablock &
|
|
|
|
, associate_blocks, deallocate_datablock
|
2011-05-05 18:25:24 -03:00
|
|
|
use blocks , only : get_nleafs
|
2011-04-21 08:39:31 -03:00
|
|
|
use error , only : print_info, print_error
|
2009-01-03 22:49:04 -06:00
|
|
|
#ifdef MPI
|
2010-12-01 09:25:30 -02:00
|
|
|
use mpitools , only : ncpus, ncpu, is_master, mallreducesuml, msendf, mrecvf
|
2009-01-03 22:49:04 -06:00
|
|
|
#endif /* MPI */
|
2010-12-01 09:25:30 -02:00
|
|
|
use problem , only : check_ref
|
2011-05-02 23:35:48 -03:00
|
|
|
use timer , only : start_timer, stop_timer
|
2010-12-01 09:25:30 -02:00
|
|
|
use variables, only : nqt
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
logical :: flag
|
2009-09-26 01:51:44 -03:00
|
|
|
integer(kind=4) :: i, j, k, l, n, p
|
2009-01-03 22:49:04 -06:00
|
|
|
|
|
|
|
#ifdef MPI
|
2011-04-23 17:29:20 -03:00
|
|
|
! tag for the MPI data exchange
|
2009-09-26 00:10:32 -03:00
|
|
|
!
|
2011-04-23 17:29:20 -03:00
|
|
|
integer(kind=4) :: itag
|
2009-09-26 00:10:32 -03:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! array for storing the refinement flags
|
2009-09-22 19:32:57 -03:00
|
|
|
!
|
2011-04-23 17:29:20 -03:00
|
|
|
integer(kind=4), dimension(:), allocatable :: ibuf
|
2009-09-26 00:10:32 -03:00
|
|
|
|
2010-12-14 11:19:11 -02:00
|
|
|
! array for number of data block for autobalancing
|
|
|
|
!
|
2011-04-23 17:29:20 -03:00
|
|
|
integer(kind=4), dimension(0:ncpus-1) :: lb
|
2010-12-14 11:19:11 -02:00
|
|
|
|
2009-09-26 00:10:32 -03:00
|
|
|
! local buffer for data block exchange
|
|
|
|
!
|
2011-04-23 17:29:20 -03:00
|
|
|
real(kind=8) , dimension(nqt,im,jm,km) :: rbuf
|
2009-01-03 22:49:04 -06:00
|
|
|
#endif /* MPI */
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_meta), pointer :: pmeta, pneigh, pchild, pparent
|
|
|
|
type(block_data), pointer :: pdata
|
|
|
|
|
2009-05-18 22:46:19 +02:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2011-05-02 23:35:48 -03:00
|
|
|
! start the mesh timer
|
|
|
|
!
|
|
|
|
call start_timer(5)
|
|
|
|
|
2011-04-21 08:39:31 -03:00
|
|
|
#ifdef DEBUG
|
|
|
|
! check mesh
|
|
|
|
!
|
|
|
|
call check_mesh('before update_mesh')
|
|
|
|
|
|
|
|
#endif /* DEBUG */
|
2011-04-23 17:29:20 -03:00
|
|
|
! iterate over elements of the data block list
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
pdata => list_data
|
|
|
|
do while (associated(pdata))
|
2011-04-23 17:29:20 -03:00
|
|
|
|
|
|
|
! assign a pointer to the meta block associated with the current data block
|
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
pmeta => pdata%meta
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! if the current data block has a meta block associated
|
|
|
|
!
|
|
|
|
if (associated(pmeta)) then
|
|
|
|
|
|
|
|
! if the associated meta block is a leaf
|
|
|
|
!
|
|
|
|
if (pmeta%leaf) then
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! check the refinement criterion for the current data block
|
|
|
|
!
|
|
|
|
pmeta%refine = check_ref(pdata)
|
2010-09-19 11:38:55 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! correct the refinement of the block for the base and top levels
|
|
|
|
!
|
|
|
|
if (pmeta%level .eq. 1) pmeta%refine = max(0, pmeta%refine)
|
|
|
|
if (pmeta%level .eq. maxlev) pmeta%refine = min(0, pmeta%refine)
|
|
|
|
|
|
|
|
end if ! pmeta is a leaf
|
|
|
|
end if ! pmeta associated
|
2009-09-21 16:57:34 -03:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! assign a pointer to the next data block
|
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
pdata => pdata%next
|
2011-04-23 17:29:20 -03:00
|
|
|
|
|
|
|
end do
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
#ifdef MPI
|
2011-04-23 17:29:20 -03:00
|
|
|
! allocate buffer for the refinement field values
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
allocate(ibuf(get_nleafs()))
|
2011-04-23 17:29:20 -03:00
|
|
|
|
|
|
|
! reset the buffer
|
|
|
|
!
|
|
|
|
ibuf(:) = 0
|
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
! store refinement flags for all blocks for exchange between processors
|
|
|
|
!
|
|
|
|
l = 1
|
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
|
|
|
if (pmeta%leaf) then
|
|
|
|
ibuf(l) = pmeta%refine
|
2011-04-23 17:29:20 -03:00
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
l = l + 1
|
|
|
|
end if
|
|
|
|
pmeta => pmeta%next
|
2009-01-08 00:08:52 -06:00
|
|
|
end do
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! update refinement flags across all processors
|
2009-09-21 16:57:34 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
call mallreducesuml(get_nleafs(), ibuf(1:get_nleafs()))
|
2009-09-21 16:57:34 -03:00
|
|
|
|
|
|
|
! update non-local block refinement flags
|
2009-09-09 14:21:46 -03:00
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
l = 1
|
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
|
|
|
if (pmeta%leaf) then
|
2011-04-23 17:29:20 -03:00
|
|
|
pmeta%refine = ibuf(l)
|
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
l = l + 1
|
|
|
|
end if
|
|
|
|
pmeta => pmeta%next
|
|
|
|
end do
|
2009-09-09 14:21:46 -03:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! deallocate the buffer
|
|
|
|
!
|
|
|
|
if (allocated(ibuf)) deallocate(ibuf)
|
|
|
|
|
|
|
|
#endif /* MPI */
|
|
|
|
! iterate over all levels starting from top and correct the refinement of blocks
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-05-18 22:46:19 +02:00
|
|
|
do l = maxlev, 1, -1
|
2009-01-08 00:08:52 -06:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! iterate over all meta blocks
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
2009-01-03 22:49:04 -06:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! check only leafs at the current level
|
|
|
|
!
|
|
|
|
if (pmeta%leaf .and. pmeta%level .eq. l) then
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! iterte over all neighbors of the current leaf
|
|
|
|
!
|
|
|
|
do i = 1, ndims
|
|
|
|
do j = 1, nsides
|
|
|
|
do k = 1, nfaces
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! assign a pointer to the current neighbor
|
|
|
|
!
|
|
|
|
pneigh => pmeta%neigh(i,j,k)%ptr
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! check if the pointer is associated with any block
|
|
|
|
!
|
|
|
|
if (associated(pneigh)) then
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
!= conditions for blocks which are selected to be refined
|
|
|
|
!
|
|
|
|
if (pmeta%refine .eq. 1) then
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! if the neighbor is set to be derefined, reset its flags (this applies to
|
|
|
|
! blocks at the current and lower levels)
|
2009-05-18 22:46:19 +02:00
|
|
|
!
|
2011-04-23 17:29:20 -03:00
|
|
|
pneigh%refine = max(0, pneigh%refine)
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! if the neighbor is at lower level, always set it to be refined
|
|
|
|
!
|
|
|
|
if (pneigh%level .lt. pmeta%level) pneigh%refine = 1
|
|
|
|
|
|
|
|
end if ! refine = 1
|
|
|
|
|
|
|
|
!= conditions for blocks which stay at the same level
|
|
|
|
!
|
|
|
|
if (pmeta%refine .eq. 0) then
|
|
|
|
|
|
|
|
! if the neighbor lays at lower level and is set to be derefined, cancel its
|
|
|
|
! derefinement
|
|
|
|
!
|
|
|
|
if (pneigh%level .lt. pmeta%level) &
|
|
|
|
pneigh%refine = max(0, pneigh%refine)
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
end if ! refine = 0
|
|
|
|
|
|
|
|
!= conditions for blocks which are selected to be derefined
|
|
|
|
!
|
|
|
|
if (pmeta%refine .eq. -1) then
|
|
|
|
|
|
|
|
! if the neighbor is at lower level and is set to be derefined, cancel its
|
|
|
|
! derefinement
|
|
|
|
!
|
|
|
|
if (pneigh%level .lt. pmeta%level) &
|
|
|
|
pneigh%refine = max(0, pneigh%refine)
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! if the neighbor is set to be refined, cancel derefinement of the current block
|
|
|
|
!
|
|
|
|
if (pneigh%refine .eq. 1) pmeta%refine = 0
|
|
|
|
|
|
|
|
end if ! refine = -1
|
|
|
|
|
|
|
|
end if ! associated(pneigh)
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
end do
|
|
|
|
end do
|
2009-01-08 00:08:52 -06:00
|
|
|
end do
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
end if ! leafs at level l
|
|
|
|
|
|
|
|
! assign a pointer to the next block
|
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
pmeta => pmeta%next
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
end do ! meta blocks
|
|
|
|
|
|
|
|
!= now check all derefined block if their siblings are set for derefinement too
|
|
|
|
! and are at the same level; check ony levels >= 2
|
2009-05-18 22:46:19 +02:00
|
|
|
!
|
2011-04-23 17:29:20 -03:00
|
|
|
if (l .ge. 2) then
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! iterate over all blocks
|
|
|
|
!
|
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! check only leafs at the current level
|
|
|
|
!
|
|
|
|
if (pmeta%leaf .and. pmeta%level .eq. l) then
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! check blocks which are selected to be derefined
|
|
|
|
!
|
|
|
|
if (pmeta%refine .eq. -1) then
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! assign a pointer to the parent of the current block
|
|
|
|
!
|
|
|
|
pparent => pmeta%parent
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! check if parent is associated with any block
|
|
|
|
!
|
|
|
|
if (associated(pparent)) then
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! reset derefinement flag
|
|
|
|
!
|
|
|
|
flag = .true.
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! iterate over all children
|
|
|
|
!
|
|
|
|
do p = 1, nchild
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2011-04-23 17:29:20 -03:00
|
|
|
! assign a pointer to the current child
|
|
|
|
!
|
|
|
|
pchild => pparent%child(p)%ptr
|
|
|
|
|
|
|
|
! check if the current child is a leaf
|
|
|
|
!
|
|
|
|
flag = flag .and. (pchild%leaf)
|
|
|
|
|
|
|
|
! check if the current child is set to be derefined
|
|
|
|
!
|
|
|
|
flag = flag .and. (pchild%refine .eq. -1)
|
|
|
|
|
|
|
|
end do ! over all children
|
|
|
|
|
|
|
|
! if not all children are proper for derefinement, cancel derefinement of all
|
|
|
|
! children
|
|
|
|
!
|
|
|
|
if (.not. flag) then
|
|
|
|
|
|
|
|
! iterate over all children
|
|
|
|
!
|
|
|
|
do p = 1, nchild
|
|
|
|
|
|
|
|
! assign a pointer to the current child
|
|
|
|
!
|
|
|
|
pchild => pparent%child(p)%ptr
|
|
|
|
|
|
|
|
! reset derefinement of the current child
|
|
|
|
!
|
|
|
|
pchild%refine = max(0, pchild%refine)
|
|
|
|
|
|
|
|
end do ! children
|
|
|
|
|
|
|
|
end if ! ~flag
|
|
|
|
|
|
|
|
end if ! pparent is associated
|
|
|
|
|
|
|
|
end if ! refine = -1
|
|
|
|
|
|
|
|
end if ! leafs at level l
|
|
|
|
|
|
|
|
! assign a pointer to the next block
|
|
|
|
!
|
|
|
|
pmeta => pmeta%next
|
|
|
|
|
|
|
|
end do ! meta blocks
|
|
|
|
|
|
|
|
end if ! l >= 2
|
|
|
|
|
|
|
|
end do ! levels
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-26 00:10:32 -03:00
|
|
|
#ifdef MPI
|
|
|
|
! find all sibling blocks which are spread over different processors
|
|
|
|
!
|
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
|
|
|
if (.not. pmeta%leaf) then
|
|
|
|
if (pmeta%child(1)%ptr%refine .eq. -1) then
|
|
|
|
|
|
|
|
! check if the parent blocks is on the same processor as the next block, if not
|
|
|
|
! move it to the same processor
|
|
|
|
!
|
2009-09-26 00:14:23 -03:00
|
|
|
if (pmeta%cpu .ne. pmeta%next%cpu) &
|
2009-09-26 00:10:32 -03:00
|
|
|
pmeta%cpu = pmeta%next%cpu
|
|
|
|
|
|
|
|
! find the case when child blocks are spread across at least 2 processors
|
|
|
|
!
|
|
|
|
flag = .false.
|
|
|
|
do p = 1, nchild
|
|
|
|
flag = flag .or. (pmeta%child(p)%ptr%cpu .ne. pmeta%cpu)
|
|
|
|
end do
|
|
|
|
|
|
|
|
if (flag) then
|
|
|
|
|
|
|
|
! iterate over all children
|
|
|
|
!
|
|
|
|
do p = 1, nchild
|
|
|
|
|
2009-09-26 01:51:44 -03:00
|
|
|
! generate the tag for communication
|
2009-09-26 00:10:32 -03:00
|
|
|
!
|
|
|
|
itag = pmeta%child(p)%ptr%cpu * ncpus + pmeta%cpu + ncpus + p + 1
|
|
|
|
|
|
|
|
! if the current children is not on the same processor, then ...
|
|
|
|
!
|
|
|
|
if (pmeta%child(p)%ptr%cpu .ne. pmeta%cpu) then
|
|
|
|
|
|
|
|
! allocate data blocks for children on the processor which will receive data
|
|
|
|
!
|
|
|
|
if (pmeta%cpu .eq. ncpu) then
|
|
|
|
call append_datablock(pdata)
|
|
|
|
call associate_blocks(pmeta%child(p)%ptr, pdata)
|
|
|
|
|
|
|
|
! receive the data
|
|
|
|
!
|
|
|
|
call mrecvf(size(rbuf), pmeta%child(p)%ptr%cpu, itag, rbuf)
|
|
|
|
|
|
|
|
! coppy buffer to data
|
|
|
|
!
|
|
|
|
pmeta%child(p)%ptr%data%u(:,:,:,:) = rbuf(:,:,:,:)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! send data to the right processor and deallocate data block
|
|
|
|
!
|
|
|
|
if (pmeta%child(p)%ptr%cpu .eq. ncpu) then
|
|
|
|
|
|
|
|
! copy data to buffer
|
|
|
|
!
|
|
|
|
rbuf(:,:,:,:) = pmeta%child(p)%ptr%data%u(:,:,:,:)
|
|
|
|
|
|
|
|
! send data
|
|
|
|
!
|
|
|
|
call msendf(size(rbuf), pmeta%cpu, itag, rbuf)
|
|
|
|
|
|
|
|
! deallocate data block
|
|
|
|
!
|
|
|
|
call deallocate_datablock(pmeta%child(p)%ptr%data)
|
2009-09-26 00:14:23 -03:00
|
|
|
end if
|
2009-09-26 00:10:32 -03:00
|
|
|
|
|
|
|
! set the current processor of the block
|
|
|
|
!
|
2009-09-26 00:14:23 -03:00
|
|
|
pmeta%child(p)%ptr%cpu = pmeta%cpu
|
2009-09-26 00:10:32 -03:00
|
|
|
end if
|
|
|
|
end do
|
|
|
|
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
pmeta => pmeta%next
|
|
|
|
end do
|
|
|
|
#endif /* MPI */
|
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
! perform the actual derefinement
|
2009-05-18 22:46:19 +02:00
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
do l = maxlev, 2, -1
|
2009-05-18 22:46:19 +02:00
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
|
|
|
|
2009-09-25 19:24:02 -03:00
|
|
|
if (pmeta%leaf) then
|
|
|
|
if (pmeta%level .eq. l) then
|
2009-09-21 16:57:34 -03:00
|
|
|
if (pmeta%refine .eq. -1) then
|
|
|
|
pparent => pmeta%parent
|
|
|
|
|
|
|
|
if (associated(pparent)) then
|
2010-02-11 23:30:46 -02:00
|
|
|
#ifdef MPI
|
2009-09-25 19:24:02 -03:00
|
|
|
if (pmeta%cpu .eq. ncpu) then
|
2010-02-11 23:30:46 -02:00
|
|
|
#endif /* MPI */
|
2009-09-25 19:24:02 -03:00
|
|
|
if (.not. associated(pparent%data)) then
|
|
|
|
call append_datablock(pdata)
|
|
|
|
call associate_blocks(pparent, pdata)
|
|
|
|
end if
|
2009-09-21 16:57:34 -03:00
|
|
|
call restrict_block(pparent)
|
2010-02-11 23:30:46 -02:00
|
|
|
#ifdef MPI
|
2009-09-25 19:24:02 -03:00
|
|
|
end if
|
2010-02-11 23:30:46 -02:00
|
|
|
#endif /* MPI */
|
2009-09-25 19:24:02 -03:00
|
|
|
|
2009-09-21 16:57:34 -03:00
|
|
|
call derefine_block(pparent)
|
|
|
|
pmeta => pparent
|
2011-04-21 08:39:31 -03:00
|
|
|
else
|
|
|
|
call print_error("mesh::update_mesh" &
|
|
|
|
, "Parent of derefined block is not associated!")
|
2009-09-22 19:32:57 -03:00
|
|
|
end if
|
2009-09-21 16:57:34 -03:00
|
|
|
end if
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
pmeta => pmeta%next
|
|
|
|
end do
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
! perform the actual refinement
|
|
|
|
!
|
|
|
|
do l = 1, maxlev - 1
|
|
|
|
|
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
|
|
|
if (pmeta%leaf) then
|
|
|
|
if (pmeta%level .eq. l) then
|
|
|
|
if (pmeta%refine .eq. 1) then
|
|
|
|
pparent => pmeta
|
2010-02-11 23:30:46 -02:00
|
|
|
#ifdef MPI
|
2009-09-21 16:57:34 -03:00
|
|
|
if (pmeta%cpu .eq. ncpu) then
|
2010-02-11 23:30:46 -02:00
|
|
|
#endif /* MPI */
|
2011-05-05 16:51:35 -03:00
|
|
|
call refine_block(pmeta, res(pmeta%level + 1,:), .true.)
|
2009-09-21 16:57:34 -03:00
|
|
|
call prolong_block(pparent)
|
2009-09-25 19:24:02 -03:00
|
|
|
call deallocate_datablock(pparent%data)
|
2010-02-11 23:30:46 -02:00
|
|
|
#ifdef MPI
|
2009-09-21 16:57:34 -03:00
|
|
|
else
|
2011-05-05 16:51:35 -03:00
|
|
|
call refine_block(pmeta, res(pmeta%level + 1,:), .false.)
|
2009-09-21 16:57:34 -03:00
|
|
|
end if
|
2010-02-11 23:30:46 -02:00
|
|
|
#endif /* MPI */
|
2009-09-21 16:57:34 -03:00
|
|
|
end if
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
pmeta => pmeta%next
|
|
|
|
end do
|
|
|
|
|
|
|
|
end do
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-22 19:32:57 -03:00
|
|
|
#ifdef MPI
|
2009-09-26 01:51:44 -03:00
|
|
|
!! AUTO BALANCING
|
2009-09-22 19:32:57 -03:00
|
|
|
!!
|
2009-09-26 01:51:44 -03:00
|
|
|
! calculate the new division
|
2009-09-22 19:32:57 -03:00
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
l = mod(get_nleafs(), ncpus) - 1
|
|
|
|
lb( : ) = get_nleafs() / ncpus
|
2010-12-14 11:19:11 -02:00
|
|
|
lb(0:l) = lb(0:l) + 1
|
2009-09-22 19:32:57 -03:00
|
|
|
|
2009-09-26 01:51:44 -03:00
|
|
|
! iterate over all metablocks and reassign the processor numbers
|
2009-09-22 19:32:57 -03:00
|
|
|
!
|
2009-09-26 01:51:44 -03:00
|
|
|
n = 0
|
2010-12-14 11:19:11 -02:00
|
|
|
l = 0
|
2009-09-22 19:32:57 -03:00
|
|
|
|
2009-09-26 01:51:44 -03:00
|
|
|
pmeta => list_meta
|
|
|
|
do while (associated(pmeta))
|
|
|
|
|
|
|
|
! assign the cpu to the current block
|
2009-09-22 19:32:57 -03:00
|
|
|
!
|
2009-09-26 01:51:44 -03:00
|
|
|
if (pmeta%cpu .ne. n) then
|
|
|
|
|
|
|
|
if (pmeta%leaf) then
|
|
|
|
|
|
|
|
! generate the tag for communication
|
|
|
|
!
|
|
|
|
itag = pmeta%cpu * ncpus + n + ncpus + 1
|
|
|
|
|
|
|
|
if (ncpu .eq. pmeta%cpu) then
|
|
|
|
|
|
|
|
! copy data to buffer
|
|
|
|
!
|
|
|
|
rbuf(:,:,:,:) = pmeta%data%u(:,:,:,:)
|
|
|
|
|
|
|
|
! send data
|
|
|
|
!
|
|
|
|
call msendf(size(rbuf), n, itag, rbuf)
|
|
|
|
|
|
|
|
! deallocate data block
|
|
|
|
!
|
|
|
|
call deallocate_datablock(pmeta%data)
|
|
|
|
|
|
|
|
! send data block
|
|
|
|
!
|
|
|
|
end if
|
|
|
|
|
|
|
|
if (ncpu .eq. n) then
|
2009-09-22 19:32:57 -03:00
|
|
|
|
2009-09-26 01:51:44 -03:00
|
|
|
! allocate data block for the current block
|
2009-09-22 19:32:57 -03:00
|
|
|
!
|
2009-09-26 01:51:44 -03:00
|
|
|
call append_datablock(pdata)
|
|
|
|
call associate_blocks(pmeta, pdata)
|
|
|
|
|
|
|
|
! receive the data
|
|
|
|
!
|
|
|
|
call mrecvf(size(rbuf), pmeta%cpu, itag, rbuf)
|
|
|
|
|
|
|
|
! coppy buffer to data
|
|
|
|
!
|
|
|
|
pmeta%data%u(:,:,:,:) = rbuf(:,:,:,:)
|
|
|
|
|
|
|
|
! receive data block
|
|
|
|
!
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
! set new processor number
|
|
|
|
!
|
|
|
|
pmeta%cpu = n
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
! increase the number of blocks on the current process; if it exceeds the
|
|
|
|
! allowed number reset the counter and increase the processor number
|
|
|
|
!
|
|
|
|
if (pmeta%leaf) then
|
2010-12-14 11:19:11 -02:00
|
|
|
l = l + 1
|
|
|
|
if (l .ge. lb(n)) then
|
2009-09-26 01:51:44 -03:00
|
|
|
n = min(ncpus - 1, n + 1)
|
2010-12-14 11:19:11 -02:00
|
|
|
l = 0
|
2009-09-26 01:51:44 -03:00
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
! assign pointer to the next block
|
|
|
|
!
|
|
|
|
pmeta => pmeta%next
|
|
|
|
end do
|
2009-09-22 19:32:57 -03:00
|
|
|
#endif /* MPI */
|
|
|
|
|
2011-04-21 08:39:31 -03:00
|
|
|
#ifdef DEBUG
|
|
|
|
! check mesh
|
|
|
|
!
|
|
|
|
call check_mesh('after update_mesh')
|
|
|
|
#endif /* DEBUG */
|
2011-05-02 23:35:48 -03:00
|
|
|
|
|
|
|
! stop the mesh timer
|
|
|
|
!
|
|
|
|
call stop_timer(5)
|
|
|
|
|
2008-12-13 15:08:18 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine update_mesh
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! prolong_block: subroutine expands the block data and copy them to children
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine prolong_block(pblock)
|
|
|
|
|
2010-12-01 09:25:30 -02:00
|
|
|
use blocks , only : block_meta, block_data, nchild
|
2010-12-05 09:46:16 -02:00
|
|
|
use config , only : ng, nh, in, jn, kn, im, jm, km
|
2010-12-05 09:56:27 -02:00
|
|
|
use interpolation, only : expand
|
2010-12-01 09:25:30 -02:00
|
|
|
use variables , only : nfl, nqt
|
|
|
|
#ifdef MHD
|
|
|
|
use variables , only : ibx, iby, ibz
|
2010-12-01 13:54:46 -02:00
|
|
|
#ifdef GLM
|
|
|
|
use variables , only : iph
|
|
|
|
#endif /* GLM */
|
2010-12-01 09:25:30 -02:00
|
|
|
#endif /* MHD */
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2009-09-29 15:32:58 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pblock
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2009-09-29 15:32:58 -03:00
|
|
|
integer :: q, p
|
|
|
|
integer :: il, iu, jl, ju, kl, ku
|
|
|
|
integer :: is, js, ks
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
2009-09-29 15:32:58 -03:00
|
|
|
integer, dimension(3) :: dm, fm
|
|
|
|
|
|
|
|
! local allocatable arrays
|
|
|
|
!
|
|
|
|
real, dimension(:,:,:,:), allocatable :: u
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_meta), pointer :: pchild
|
2010-09-19 11:46:00 +02:00
|
|
|
type(block_data), pointer :: pdata
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2010-09-19 11:46:00 +02:00
|
|
|
! assign the pdata pointer
|
|
|
|
!
|
|
|
|
pdata => pblock%data
|
|
|
|
|
2009-09-29 15:32:58 -03:00
|
|
|
! prepare dimensions
|
|
|
|
!
|
|
|
|
dm(:) = (/ im, jm, km /)
|
2008-12-13 15:08:18 -06:00
|
|
|
fm(:) = 2 * (dm(:) - ng)
|
2009-09-29 15:32:58 -03:00
|
|
|
#if NDIMS == 2
|
|
|
|
fm(3) = 1
|
|
|
|
ks = 1
|
|
|
|
kl = 1
|
|
|
|
ku = 1
|
|
|
|
#endif /* NDIMS == 2 */
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-29 15:32:58 -03:00
|
|
|
! allocate array to the product of expansion
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2010-02-28 18:35:57 -03:00
|
|
|
allocate(u(nqt, fm(1), fm(2), fm(3)))
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-29 15:32:58 -03:00
|
|
|
! expand all variables and place them in the array u
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2010-03-08 19:09:49 -03:00
|
|
|
do q = 1, nfl
|
2010-12-05 09:56:27 -02:00
|
|
|
call expand(dm(:), fm(:), nh, pdata%u(q,:,:,:), u(q,:,:,:))
|
2008-12-13 15:08:18 -06:00
|
|
|
end do
|
|
|
|
|
2010-03-08 19:09:49 -03:00
|
|
|
#ifdef MHD
|
2010-12-01 13:54:46 -02:00
|
|
|
! expand the cell centered magnetic field components
|
|
|
|
!
|
|
|
|
do q = ibx, ibz
|
2010-12-05 09:56:27 -02:00
|
|
|
call expand(dm(:), fm(:), nh, pdata%u(q,:,:,:), u(q,:,:,:))
|
2010-12-01 13:54:46 -02:00
|
|
|
end do
|
2010-12-05 09:56:27 -02:00
|
|
|
#ifdef GLM
|
2010-12-01 13:54:46 -02:00
|
|
|
|
|
|
|
! expand the scalar potential Psi
|
|
|
|
!
|
2010-12-05 09:56:27 -02:00
|
|
|
call expand(dm(:), fm(:), nh, pdata%u(iph,:,:,:), u(iph,:,:,:))
|
2010-12-01 13:54:46 -02:00
|
|
|
#endif /* GLM */
|
2010-03-08 19:09:49 -03:00
|
|
|
#endif /* MHD */
|
2009-09-29 15:32:58 -03:00
|
|
|
! iterate over all children
|
|
|
|
!
|
|
|
|
do p = 1, nchild
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-29 15:32:58 -03:00
|
|
|
! assign pointer to the current child
|
|
|
|
!
|
|
|
|
pchild => pblock%child(p)%ptr
|
|
|
|
|
2010-12-05 10:14:57 -02:00
|
|
|
! obtain the position of child in the parent block
|
2009-09-29 15:32:58 -03:00
|
|
|
!
|
2010-12-05 10:14:57 -02:00
|
|
|
is = pchild%pos(1)
|
|
|
|
js = pchild%pos(2)
|
2009-09-29 15:32:58 -03:00
|
|
|
#if NDIMS == 3
|
2010-12-05 10:14:57 -02:00
|
|
|
ks = pchild%pos(3)
|
2009-09-29 15:32:58 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
! calculate indices of the current child subdomain
|
|
|
|
!
|
|
|
|
il = 1 + is * in
|
|
|
|
jl = 1 + js * jn
|
|
|
|
#if NDIMS == 3
|
|
|
|
kl = 1 + ks * kn
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
iu = il + im - 1
|
|
|
|
ju = jl + jm - 1
|
|
|
|
#if NDIMS == 3
|
|
|
|
ku = kl + km - 1
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
|
|
|
|
! copy data to the current child
|
|
|
|
!
|
2010-02-28 18:35:57 -03:00
|
|
|
pchild%data%u(1:nqt,1:im,1:jm,1:km) = u(1:nqt,il:iu,jl:ju,kl:ku)
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
end do
|
|
|
|
|
2009-09-29 15:32:58 -03:00
|
|
|
! deallocate local arrays
|
|
|
|
!
|
|
|
|
if (allocated(u)) deallocate(u)
|
|
|
|
|
2008-12-13 15:08:18 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine prolong_block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! restrict_block: subroutine shrinks the block data and copy them from children
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine restrict_block(pblock)
|
|
|
|
|
2010-12-01 09:25:30 -02:00
|
|
|
use blocks , only : block_meta, block_data, nchild
|
2010-09-19 11:49:21 +02:00
|
|
|
use config , only : ng, in, ih, im, ib, ie, nh, jn, jh, jm, jb, je &
|
|
|
|
, kn, kh, km, kb, ke
|
2010-12-01 09:25:30 -02:00
|
|
|
use variables , only : nfl
|
|
|
|
#ifdef MHD
|
|
|
|
use variables , only : ibx, iby, ibz
|
2010-12-01 13:45:54 -02:00
|
|
|
#ifdef GLM
|
|
|
|
use variables , only : iph
|
|
|
|
#endif /* GLM */
|
2010-12-01 09:25:30 -02:00
|
|
|
#endif /* MHD */
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2009-09-21 16:57:34 -03:00
|
|
|
type(block_meta), pointer, intent(inout) :: pblock
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2010-09-19 11:49:21 +02:00
|
|
|
integer :: p
|
|
|
|
integer :: if, jf, kf
|
|
|
|
integer :: il, jl, kl, iu, ju, ku
|
|
|
|
integer :: ip, jp, kp
|
|
|
|
integer :: is, js, ks, it, jt, kt
|
2010-02-11 23:30:46 -02:00
|
|
|
|
2008-12-13 15:08:18 -06:00
|
|
|
! local pointers
|
|
|
|
!
|
2010-09-19 11:49:21 +02:00
|
|
|
type(block_data), pointer :: pparent, pchild
|
2008-12-13 15:08:18 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2010-09-19 11:49:21 +02:00
|
|
|
! assign pointers
|
2010-02-11 23:30:46 -02:00
|
|
|
!
|
2010-09-19 11:49:21 +02:00
|
|
|
pparent => pblock%data
|
2010-02-11 23:30:46 -02:00
|
|
|
|
2009-09-29 15:12:19 -03:00
|
|
|
! iterate over all children
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2009-09-29 15:12:19 -03:00
|
|
|
do p = 1, nchild
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2009-09-29 15:12:19 -03:00
|
|
|
! assign pointer to the current child
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2010-09-19 11:49:21 +02:00
|
|
|
pchild => pblock%child(p)%ptr%data
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2010-12-05 10:16:17 -02:00
|
|
|
! obtain the position of the current child in the parent block
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2010-12-05 10:19:32 -02:00
|
|
|
if = pchild%meta%pos(1)
|
|
|
|
jf = pchild%meta%pos(2)
|
2009-09-29 15:12:19 -03:00
|
|
|
#if NDIMS == 3
|
2010-12-05 10:19:32 -02:00
|
|
|
kf = pchild%meta%pos(3)
|
2009-09-29 15:12:19 -03:00
|
|
|
#endif /* NDIMS == 3 */
|
2008-12-13 15:08:18 -06:00
|
|
|
|
2010-09-19 11:49:21 +02:00
|
|
|
! calculate the bound indices of the source nad destination arrays
|
|
|
|
!
|
|
|
|
if (if .eq. 0) then
|
|
|
|
il = 1
|
|
|
|
iu = ie
|
|
|
|
is = ib - nh
|
|
|
|
it = ih
|
|
|
|
else
|
|
|
|
il = ib
|
|
|
|
iu = im
|
|
|
|
is = ih + 1
|
|
|
|
it = ie + nh
|
|
|
|
end if
|
|
|
|
ip = il + 1
|
|
|
|
if (jf .eq. 0) then
|
|
|
|
jl = 1
|
|
|
|
ju = je
|
|
|
|
js = jb - nh
|
|
|
|
jt = jh
|
|
|
|
else
|
|
|
|
jl = jb
|
|
|
|
ju = jm
|
|
|
|
js = jh + 1
|
|
|
|
jt = je + nh
|
|
|
|
end if
|
|
|
|
jp = jl + 1
|
|
|
|
#if NDIMS == 3
|
|
|
|
if (kf .eq. 0) then
|
|
|
|
kl = 1
|
|
|
|
ku = ke
|
|
|
|
ks = kb - nh
|
|
|
|
kt = kh
|
|
|
|
else
|
|
|
|
kl = kb
|
|
|
|
ku = km
|
|
|
|
ks = kh + 1
|
|
|
|
kt = ke + nh
|
|
|
|
end if
|
|
|
|
kp = kl + 1
|
|
|
|
#endif /* NDIMS == 3 */
|
2010-02-11 23:30:46 -02:00
|
|
|
|
2010-09-19 11:49:21 +02:00
|
|
|
! copy the variables from the current child to the proper location of
|
|
|
|
! the parent block
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
2010-07-04 01:15:34 -03:00
|
|
|
#if NDIMS == 2
|
2010-09-19 11:49:21 +02:00
|
|
|
pparent%u(1:nfl,is:it,js:jt,1) = &
|
|
|
|
0.25 * (pchild%u(1:nfl,il:iu:2,jl:ju:2,1) &
|
|
|
|
+ pchild%u(1:nfl,ip:iu:2,jl:ju:2,1) &
|
|
|
|
+ pchild%u(1:nfl,il:iu:2,jp:ju:2,1) &
|
|
|
|
+ pchild%u(1:nfl,ip:iu:2,jp:ju:2,1))
|
2010-02-11 23:30:46 -02:00
|
|
|
|
|
|
|
#ifdef MHD
|
2010-12-01 13:45:54 -02:00
|
|
|
pparent%u(ibx:ibz,is:it,js:jt,1) = &
|
|
|
|
0.25 * (pchild%u(ibx:ibz,il:iu:2,jl:ju:2,1) &
|
|
|
|
+ pchild%u(ibx:ibz,ip:iu:2,jl:ju:2,1) &
|
|
|
|
+ pchild%u(ibx:ibz,il:iu:2,jp:ju:2,1) &
|
|
|
|
+ pchild%u(ibx:ibz,ip:iu:2,jp:ju:2,1))
|
2010-12-14 16:34:51 -02:00
|
|
|
#ifdef GLM
|
2010-12-01 13:45:54 -02:00
|
|
|
pparent%u(iph ,is:it,js:jt,1) = &
|
|
|
|
0.25 * (pchild%u(iph ,il:iu:2,jl:ju:2,1) &
|
|
|
|
+ pchild%u(iph ,ip:iu:2,jl:ju:2,1) &
|
|
|
|
+ pchild%u(iph ,il:iu:2,jp:ju:2,1) &
|
|
|
|
+ pchild%u(iph ,ip:iu:2,jp:ju:2,1))
|
|
|
|
#endif /* GLM */
|
2010-02-11 23:30:46 -02:00
|
|
|
#endif /* MHD */
|
2010-07-04 01:15:34 -03:00
|
|
|
#endif /* NDIMS == 2 */
|
2010-09-19 11:49:21 +02:00
|
|
|
#if NDIMS == 3
|
|
|
|
pparent%u(1:nfl,is:it,js:jt,ks:kt) = &
|
|
|
|
0.125 * (pchild%u(1:nfl,il:iu:2,jl:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(1:nfl,ip:iu:2,jl:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(1:nfl,il:iu:2,jp:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(1:nfl,ip:iu:2,jp:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(1:nfl,il:iu:2,jl:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(1:nfl,ip:iu:2,jl:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(1:nfl,il:iu:2,jp:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(1:nfl,ip:iu:2,jp:ju:2,kp:ku:2))
|
|
|
|
#ifdef MHD
|
2010-12-01 13:45:54 -02:00
|
|
|
pparent%u(ibx:ibz,is:it,js:jt,ks:kt) = &
|
|
|
|
0.125 * (pchild%u(ibx:ibz,il:iu:2,jl:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(ibx:ibz,ip:iu:2,jl:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(ibx:ibz,il:iu:2,jp:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(ibx:ibz,ip:iu:2,jp:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(ibx:ibz,il:iu:2,jl:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(ibx:ibz,ip:iu:2,jl:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(ibx:ibz,il:iu:2,jp:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(ibx:ibz,ip:iu:2,jp:ju:2,kp:ku:2))
|
2010-12-14 16:34:51 -02:00
|
|
|
#ifdef GLM
|
2010-12-01 13:45:54 -02:00
|
|
|
pparent%u(iph ,is:it,js:jt,ks:kt) = &
|
|
|
|
0.125 * (pchild%u(iph ,il:iu:2,jl:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(iph ,ip:iu:2,jl:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(iph ,il:iu:2,jp:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(iph ,ip:iu:2,jp:ju:2,kl:ku:2) &
|
|
|
|
+ pchild%u(iph ,il:iu:2,jl:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(iph ,ip:iu:2,jl:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(iph ,il:iu:2,jp:ju:2,kp:ku:2) &
|
|
|
|
+ pchild%u(iph ,ip:iu:2,jp:ju:2,kp:ku:2))
|
|
|
|
#endif /* GLM */
|
2010-09-19 11:49:21 +02:00
|
|
|
#endif /* MHD */
|
|
|
|
#endif /* NDIMS == 3 */
|
|
|
|
end do
|
2008-12-13 15:08:18 -06:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine restrict_block
|
|
|
|
!
|
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
|
|
|
! clears_mesh: subroutine deallocates mesh, removing blocks
|
|
|
|
!
|
2008-12-13 15:08:18 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2011-02-27 13:49:26 -03:00
|
|
|
subroutine clear_mesh()
|
2008-11-11 16:12:26 -06:00
|
|
|
|
2011-04-28 12:04:22 -03:00
|
|
|
use blocks , only : clear_blocks
|
|
|
|
use error , only : print_info
|
|
|
|
use mpitools, only : is_master
|
2011-05-02 23:46:16 -03:00
|
|
|
use timer , only : start_timer, stop_timer
|
2008-11-11 16:12:26 -06:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2008-12-13 15:08:18 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
2011-05-02 23:46:16 -03:00
|
|
|
! start the mesh timer
|
|
|
|
!
|
|
|
|
call start_timer(5)
|
|
|
|
|
2008-11-11 16:12:26 -06:00
|
|
|
! deallocate block structure
|
|
|
|
!
|
|
|
|
call clear_blocks
|
|
|
|
|
2008-12-08 20:03:01 -06:00
|
|
|
! deallocating coordinate variables
|
|
|
|
!
|
2011-03-02 15:01:08 -03:00
|
|
|
if (allocated(ax) ) deallocate(ax)
|
|
|
|
if (allocated(ay) ) deallocate(ay)
|
|
|
|
if (allocated(az) ) deallocate(az)
|
|
|
|
if (allocated(adx) ) deallocate(adx)
|
|
|
|
if (allocated(ady) ) deallocate(ady)
|
|
|
|
if (allocated(adz) ) deallocate(adz)
|
|
|
|
if (allocated(adxi) ) deallocate(adxi)
|
|
|
|
if (allocated(adyi) ) deallocate(adyi)
|
|
|
|
if (allocated(adzi) ) deallocate(adzi)
|
|
|
|
if (allocated(advol)) deallocate(advol)
|
2011-05-05 16:51:35 -03:00
|
|
|
if (allocated(res) ) deallocate(res)
|
2008-12-08 20:03:01 -06:00
|
|
|
|
2011-04-28 12:04:22 -03:00
|
|
|
! close the handler of the mesh statistics file
|
|
|
|
!
|
|
|
|
if (is_master()) close(funit)
|
|
|
|
|
2011-05-02 23:46:16 -03:00
|
|
|
! stop the mesh timer
|
|
|
|
!
|
|
|
|
call stop_timer(5)
|
|
|
|
|
2008-12-13 15:08:18 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
|
|
|
end subroutine clear_mesh
|
2011-04-28 12:04:22 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! store_mesh_stats: subroutine stores mesh statistics
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine store_mesh_stats(n, t)
|
|
|
|
|
|
|
|
use blocks , only : block_meta, list_meta
|
2011-05-05 18:25:24 -03:00
|
|
|
use blocks , only : get_mblocks, get_nleafs
|
2011-04-28 12:04:22 -03:00
|
|
|
use config , only : ncells, nghost, maxlev
|
|
|
|
use mpitools, only : is_master, ncpus
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! arguments
|
|
|
|
!
|
|
|
|
integer , intent(in) :: n
|
|
|
|
real(kind=8), intent(in) :: t
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer(kind=4) :: l
|
|
|
|
real(kind=4) :: cov, eff
|
|
|
|
|
|
|
|
! local arrays
|
|
|
|
!
|
|
|
|
integer(kind=4), dimension(maxlev) :: ldist
|
|
|
|
#ifdef MPI
|
|
|
|
integer(kind=4), dimension(ncpus) :: cdist
|
|
|
|
#endif /* MPI */
|
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
|
|
|
type(block_meta), pointer :: pmeta
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! store the statistics about mesh
|
|
|
|
!
|
|
|
|
if (is_master()) then
|
|
|
|
|
|
|
|
! calculate the coverage
|
|
|
|
!
|
2011-05-05 18:25:24 -03:00
|
|
|
cov = (1.0 * get_nleafs()) / tblocks
|
|
|
|
eff = (1.0 * get_nleafs() * (ncells + 2 * nghost)**NDIMS) &
|
2011-04-28 12:04:22 -03:00
|
|
|
/ product(effres(1:NDIMS) + 2 * nghost)
|
|
|
|
|
|
|
|
! get the block level distribution
|
|
|
|
!
|
|
|
|
ldist(:) = 0
|
|
|
|
#ifdef MPI
|
|
|
|
cdist(:) = 0
|
|
|
|
#endif /* MPI */
|
|
|
|
pmeta => list_meta
|
|
|
|
do while(associated(pmeta))
|
|
|
|
if (pmeta%leaf) then
|
|
|
|
ldist(pmeta%level) = ldist(pmeta%level) + 1
|
|
|
|
#ifdef MPI
|
|
|
|
cdist(pmeta%cpu+1) = cdist(pmeta%cpu+1) + 1
|
|
|
|
#endif /* MPI */
|
|
|
|
end if
|
|
|
|
pmeta => pmeta%next
|
|
|
|
end do
|
|
|
|
|
|
|
|
! write down the statistics
|
|
|
|
!
|
|
|
|
write(funit,"(2x,i8,2x,1pe14.8,2(2x,i6),2(2x,1pe14.8),$)") &
|
2011-05-05 18:25:24 -03:00
|
|
|
n, t, get_nleafs(), get_mblocks(), cov, eff
|
2011-04-28 12:04:22 -03:00
|
|
|
write(funit,"(' ',$)")
|
|
|
|
do l = 1, maxlev
|
|
|
|
write(funit,"(2x,i6,$)") ldist(l)
|
|
|
|
end do
|
|
|
|
#ifdef MPI
|
|
|
|
write(funit,"(' ',$)")
|
|
|
|
do l = 1, ncpus
|
|
|
|
write(funit,"(2x,i6,$)") cdist(l)
|
|
|
|
end do
|
|
|
|
#endif /* MPI */
|
|
|
|
write(funit,"('')")
|
|
|
|
|
|
|
|
end if
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine store_mesh_stats
|
2011-04-21 08:39:31 -03:00
|
|
|
#ifdef DEBUG
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! check_mesh: subroutine checks if the block structure is correct
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2011-04-21 14:13:25 -03:00
|
|
|
subroutine check_mesh(string)
|
2011-04-21 08:39:31 -03:00
|
|
|
|
|
|
|
use blocks, only : block_meta, list_meta
|
|
|
|
use blocks, only : last_id, nchild, ndims, nsides, nfaces
|
2011-04-21 14:13:25 -03:00
|
|
|
use blocks, only : check_metablock
|
2011-04-21 08:39:31 -03:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2011-04-21 14:13:25 -03:00
|
|
|
character(len=*), intent(in) :: string
|
2011-04-21 08:39:31 -03:00
|
|
|
|
|
|
|
! local pointers
|
|
|
|
!
|
2011-04-21 14:13:25 -03:00
|
|
|
type(block_meta), pointer :: pmeta
|
2011-04-21 08:39:31 -03:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! check meta blocks
|
|
|
|
!
|
|
|
|
pmeta => list_meta
|
|
|
|
do while(associated(pmeta))
|
|
|
|
|
2011-04-21 14:13:25 -03:00
|
|
|
! check the current block
|
2011-04-21 08:39:31 -03:00
|
|
|
!
|
2011-04-21 14:13:25 -03:00
|
|
|
call check_metablock(pmeta, string)
|
2011-04-21 08:39:31 -03:00
|
|
|
|
|
|
|
pmeta => pmeta%next
|
|
|
|
end do
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine check_mesh
|
|
|
|
#endif /* DEBUG */
|
2008-11-11 16:12:26 -06:00
|
|
|
|
2008-12-13 15:08:18 -06:00
|
|
|
!===============================================================================
|
2008-11-11 16:12:26 -06:00
|
|
|
!
|
|
|
|
end module
|