[FLASH-USERS] Using guard cell corners... bad?
Aaron Jackson
Aaron.Jackson at stonybrook.edu
Thu Mar 3 13:40:02 EST 2011
Hi Kevin,
Sorry, I should clarify my last statement about the data being
incorrect. That is just my suspicion inferred from not having this
problem with a uniform grid and that FLASH blocks only seem to know
about their edge-adjacent neighbors, at least according to the original
paper. After digging around a bit more, it seems like in
amr_1blk_guardcell_srl.F90 that there is code to communicate with
diagonal-adjacent neighbors since it pulls the remote pe from a matrix
surrblks. Although, what actually happens from there is hard to
determine from looking at the code.
if(ndim.eq.3) then
!
! Loop over the 8 corners
do kk = 1,3,2
do jj = 1,3,2
do ii = 1,3,2
jpolar = 0
remote_blk = surrblks(1,ii,jj,kk)
remote_pe = surrblks(2,ii,jj,kk)
remote_type = surrblks(3,ii,jj,kk)
-Aaron
Aaron Jackson wrote:
> Hi Kevin,
>
> Thanks for the help... Currently, though diagonals are set to true, as
> far as I can tell. In the amr_runtime_parameters.dump file from my
> simulation, diagonals is set to true. The issue isn't that the corner
> cells aren't being filled, it just seems that they're being filled
> from extrapolating information in the edge-adjacent neighbors rather
> than being copied/prolonged/restricted with data from the
> diagonally-adjacent neighbor. So the data is there, it's just not
> correct.
>
> Grep for diagonals returns:
> [ajackson at ctsv obj_Flame3StageNoise]$ grep 'diagonals' *.F90
> amr_set_runtime_parameters.F90: read (35,*) diagonals
> gr_amr_dump_runtime_parameters.F90: write (35,*)
> diagonals ,', diagonals'
> mpi_amr_1blk_restrict.F90: if (.not.diagonals) then
> mpi_amr_1blk_restrict.F90: write(*,*) 'amr_1blk_restrict:
> diagonals off'
> mpi_amr_1blk_restrict.F90: ldiag = diagonals
> mpi_amr_guardcell.F90: If (.not.diagonals) then
> mpi_amr_guardcell.F90: write(*,*) 'amr_guardcell: diagonals off'
> mpi_amr_guardcell.F90: ldiag = diagonals
> mpi_amr_prolong.F90: ldiag = diagonals
> mpi_amr_restrict_fulltree.F90: if (.not.diagonals) then
> mpi_amr_restrict_fulltree.F90: write(*,*) 'amr_1blk_restrict:
> diagonals off'
> mpi_amr_restrict_fulltree.F90: ldiag = diagonals
> physicaldata.F90: public :: diagonals
> physicaldata.F90: logical, save :: diagonals
> physicaldata.F90: logical, parameter :: diagonals = .true.
> physicaldata.F90: logical, parameter :: diagonals = .false.
> test_multigrid.F90: if (diagonals) then
> test_multigrid.F90: write(*,*) 'diagonals on '
> test_multigrid.F90: ldiag = diagonals
>
>
> A piece of code from mpi_amr_guardcells.F90
>
> ldiag = diagonals
> l_srl_only = .false. ! fill srl and coarse
> icoord = 0 ! fill in all coord
> directions
> Call amr_1blk_guardcell(mype,iopt,nlayers,lb,mype, &
> lcc,lfc,lec,lnc, &
> l_srl_only,icoord,ldiag, &
> nlayersx,nlayersy,nlayersz)
>
> Thanks,
> -Aaron
>
> Kevin Olson wrote:
>> Dear Aaron,
>>
>> There is a way to force PARAMESH to fill the corners cells correctly.
>> Here is how it is done in PARAMESH:
>>
>> mype = local process no.
>> iopt = 1
>> nguard = total no. of guardcells
>>
>> diagonals = .true. ! This is what tells the guardcell filling to
>> ensure that diagonal cells are filled correctly
>> call amr_guardcell (mype, iopt, nguard)
>>
>> The key here is to set the 'diagonals' variables to be true.
>>
>> You may need to dig into the FLASH code a bit to find exactly where
>> and how they ulitmately call the guardcell filling routine and modify
>> that.
>> Note that this will not work with PARAMESH v2 !!!
>>
>> Best,
>> Kevin Olson
>>
>> On Mar 3, 2011, at 10:40 AM, Aaron Jackson wrote:
>>
>>
>>> Dear flash-users,
>>>
>>> I'm implementing a turbulence operator described originally in Colin
>>> et al. (2000) in which I need to perform a curl of the laplacian of
>>> the velocity field. This requires two finite difference operators,
>>> which to fit onto one block, requires the use of the guard cell
>>> corners. In quiescent test cases I have run, significant numerical
>>> noise develops along block boundaries, even at the same refinement
>>> level. Since I don't have this problem using a uniform grid, my
>>> guess is that corner guard cell information is not actually copied
>>> from the diagonally adjacent neighbor, but rather interpolated from
>>> the edge-adjacent neighbors.
>>>
>>> My main question is:
>>> Is there any way to get corner guard-cells filled by
>>> copying/prolonging/restricting data from the diagonally-adjacent
>>> neighbor?
>>>
>>> I know that the obvious fix is to perform the laplacian operator
>>> first, then perform a guard cell fill and then perform the curl, but
>>> in order for my turbulence measure to be meaningful, it must be
>>> associated with a particular length scale. This becomes a problem if
>>> the edge-adjacent neighbor is at a different refinement level. I've
>>> thought of introducing a different stride in the finite differences
>>> (possible odd-even decoupling), or potentially reducing the order of
>>> the finite-difference (5 to 3-point stencil) if the refinement level
>>> of the neighbor is 1 lower than maximum. But I'm uncertain about the
>>> accuracy of these choices. I would much prefer data directly from
>>> the diagonally adjacent neighbors to be filled to the guard cell
>>> corners.
>>>
>>> Thanks,
>>> -Aaron
>>>
>>> --
>>> Aaron P. Jackson
>>> Department of Physics and Astronomy
>>> Stony Brook University
>>> Stony Brook, NY 11789-3800
>>>
>>> email: Aaron.Jackson at stonybrook.edu
>>> web: http://www.astro.sunysb.edu/ajackson
>>>
>>>
>>
>>
>
--
Aaron P. Jackson
Department of Physics and Astronomy
Stony Brook University
Stony Brook, NY 11789-3800
email: Aaron.Jackson at stonybrook.edu
web: http://www.astro.sunysb.edu/ajackson
More information about the flash-users
mailing list