[FLASH-USERS] Multigrid boundaries with PFFT top level solver

Joshua Wall joshua.e.wall at gmail.com
Fri Jun 17 15:47:06 EDT 2016


Dear users

    I've been attempting to do a mixed boundary simulation that has
different numbers of root grids in each direction. I think this can be
accomplished with the multigrid solver combined with the top level PFFT
solver (Hybrid I think its called in the user guide). But I've noticed in
Gravity_getPotentialListOfBlocks for the MG solver that it sets all the
boundary types equal to the type which grav_boundary_type is set to in the
par file:

  if(.not.useGravity) return

  if(.not.updateGravity) return

  call Timers_start("gravity Barrier")
  call MPI_Barrier (grv_meshComm, ierr)
  call Timers_stop("gravity Barrier")

  call Timers_start("gravity")

  bcTypes = grav_boundary
  where (bcTypes == PERIODIC)
     bcTypes = GRID_PDE_BND_PERIODIC
  elsewhere (bcTypes == ISOLATED)
     bcTypes = GRID_PDE_BND_ISOLATED
  elsewhere (bcTypes == DIRICHLET)
     bcTypes = GRID_PDE_BND_DIRICHLET
  elsewhere (bcTypes == OUTFLOW)
     bcTypes = GRID_PDE_BND_NEUMANN
  end where
  bcValues = 0.



Then later on in hgSolve it compares the boundary type thats passed with
the type initialized from hgPFFTinit.F90 (which reads xl_mg_boundary_type,
xr, yl, yr, etc). When it sees these are different, it throws an error:


  do eachBoundary = 1, 2*NDIM
!!$     if (bndTypes(eachBoundary) == gr_hgbcTypes(eachBoundary) .OR.
suppressPfft) then
!!$        gr_hgBndTypes(eachBoundary) = bndTypes(eachBoundary)
     if (bndTypes(eachBoundary) == gr_hgbcTypes(eachBoundary) .OR.
suppressPfft .OR. &
          (bndTypes(eachBoundary)==MG_BND_GIVENVAL .AND.
gr_hgbcTypes(eachBoundary)==MG_BND_DIRICHLET) ) then
        gr_hgBndTypes(eachBoundary) = bndTypes(eachBoundary)
     else
        if (gr_meshMe .eq. 0) then
           write(*,*) 'gr_hgSolve Error: Boundary Conditions for Poisson
Solver is inconsistent.'
           write(*,*) 'gr_hgSolve Error: direction=',eachBoundary
           write(*,*) 'gr_hgSolve Error: gr_hgbcTypes(direction)
=',gr_hgbcTypes(eachBoundary)
           write(*,*) 'gr_hgSolve Error: bndTypes(direction)
=',bndTypes(eachBoundary)
        endif
        call Driver_abortFlash('gr_hgSolve Error: BC type argument
inconsistent')
     end if
  end do

So my question is can I safely change the setting in
Gravity_getPotentialListOfBlocks to read the six different boundaries used
by hgSolve and read by hgPFFTinit, since it appears it should be alright to
solve using different ones (and its mentioned in the user guide that it can
do mixed boundaries)? I ask because it seems almost as if someone changed
it to only do one type of boundary because there may have been an error
before, and I'm always a bit cautious to change such a integral piece of
code as this.

For clarity, here's my par file settings for the domain and gravity for
what I'd like to do:


 27 xmin        = -1.29093e+17 # effective resolution 32^3
 28 xmax        =  1.29093e+17
 29 ymin        = -1.29093e+17
 30 ymax        =  1.29093e+17
 31 zmin        = -1.9364e17 #-1.29093e+17
 32 zmax        =  1.9364e17 #1.29093e+17
 33
 34 nblockx = 1
 35 nblocky = 1
 36 nblockz = 2
 37
 38 geometry = "cartesian"
 39 xl_boundary_type  = "periodic" # "reflecting"
 40 xr_boundary_type  = "periodic" # "reflecting"
 41 yl_boundary_type  = "periodic" # "reflecting"
 42 yr_boundary_type  = "periodic" # "reflecting"
 43 zl_boundary_type  = "outflow" #"periodic" # "reflecting"
 44 zr_boundary_type  = "outflow" #"periodic" # "reflecting"

...

141 useGravity         = .true.
142 grav_boundary_type = "periodic" # "isolated"
143 mg_maxResidualNorm = 1.0e-6
144 mg_printNorm       = .false.
145
146 xl_mg_boundary_type  = "periodic"
147 xr_mg_boundary_type  = "periodic"
148 yl_mg_boundary_type  = "periodic"
149 yr_mg_boundary_type  = "periodic"
150 zl_mg_boundary_type  = "outflow" #"periodic"
151 zr_mg_boundary_type  = "outflow" #"periodic"

Thanks for the assistance,

Joshua Wall
-- 
Joshua Wall
Doctoral Candidate
Department of Physics
Drexel University
3141 Chestnut Street
Philadelphia, PA 19104
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://flash.rochester.edu/pipermail/flash-users/attachments/20160617/01334809/attachment.htm>


More information about the flash-users mailing list