Skip to content

Commit

Permalink
(phantom2pdf) allocate memory for AMR grid dynamically
Browse files Browse the repository at this point in the history
  • Loading branch information
danieljprice committed Dec 9, 2019
1 parent f8eac1f commit 7c8d740
Showing 1 changed file with 19 additions and 5 deletions.
24 changes: 19 additions & 5 deletions src/utils/adaptivemesh.f90
Original file line number Diff line number Diff line change
Expand Up @@ -23,7 +23,7 @@
!+
!--------------------------------------------------------------------------
module adaptivemesh
use dim, only:maxp,maxp_hard,periodic
use dim, only:periodic
implicit none
!--this controls the number of cells on each level
integer, parameter :: nsub = 2
Expand All @@ -36,7 +36,7 @@ module adaptivemesh
! bear in mind that total number of cells is maxmeshes*nsub**ndim
! so use maxp/nsub**ndim to get number of cells = maxp
!
integer, parameter :: maxmeshes = 5e6 ! hardwired, bad but emit warning if too small
integer :: maxmeshes = 0 ! zero until allocated
!
!--resolution of the root grid (2^ifirstlevel)^ndim
!
Expand All @@ -47,15 +47,15 @@ module adaptivemesh
!--grid array stores children in first 8 cells, then parent
! (currently we do not store the parent as it is not needed)
!
integer :: gridnodes(maxchildren,maxmeshes)
integer, allocatable :: gridnodes(:,:)
! integer, parameter :: ilocparent = maxchildren + 1

integer :: maxlevel

contains

subroutine build_mesh(xyzh,np,nmesh,xmin,dxmax)
real, intent(in) :: xyzh(4,maxp)
real, intent(in) :: xyzh(:,:)
integer, intent(in) :: np
integer, intent(out) :: nmesh
real, intent(in), optional :: xmin(ndim),dxmax(ndim)
Expand Down Expand Up @@ -89,15 +89,23 @@ subroutine build_mesh(xyzh,np,nmesh,xmin,dxmax)
!
maxlevel = ifirstlevel
nmesh = 1
!
!--allocate memory
!
maxmeshes = 2*np + 1 ! allow twice the number of particles
allocate(gridnodes(maxchildren,maxmeshes))
gridnodes(:,1) = -1

!
!--build the mesh
!
!$omp parallel do default(none) schedule(guided,10) &
!$omp shared(np,xyzh,xminp,dxmaxp,nmesh) &
!$omp private(i)
do i=1,np
call refine_mesh(xyzh(:,i),1,ifirstlevel,xminp,dxmaxp,nmesh)
enddo
!$omp end parallel do

print "(a,i10,a,i10,a,i2,a,i6,a)",&
' build_mesh: nmeshes = ',nmesh,', ncells = ',nmesh*(nsub**ndim),', max level = ',maxlevel, &
' (effective resolution = ',nsub**maxlevel,'^3)'
Expand Down Expand Up @@ -242,4 +250,10 @@ recursive subroutine refine_mesh(xyzhi,imesh,level,xminl,dxmax,nmesh)

end subroutine refine_mesh

subroutine free_mesh()

if (allocated(gridnodes)) deallocate(gridnodes)

end subroutine free_mesh

end module adaptivemesh

0 comments on commit 7c8d740

Please sign in to comment.