COORDINATES: Use ncells instead of in, jn, and kn.

Also rename ng to nghosts.

Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
This commit is contained in:
Grzegorz Kowal 2019-02-04 14:43:10 -02:00
parent d4059b8075
commit af2c7213f3
6 changed files with 163 additions and 183 deletions

View File

@ -555,7 +555,7 @@ module boundaries
#endif /* MPI */
use blocks , only : ndims, nsides
use coordinates , only : minlev, maxlev
use coordinates , only : in, jn, kn
use coordinates , only : ni => ncells
use coordinates , only : ib, ie, ibl
use coordinates , only : jb, je, jbl
use coordinates , only : kb, ke, kbl
@ -607,15 +607,15 @@ module boundaries
! calculate half sizes
!
ih = in / 2
jh = jn / 2
ih = ni / 2
jh = ni / 2
#if NDIMS == 2
kh = kn
kh = 1
kl = 1
ku = kn
ku = 1
#endif /* NDIMS == 2 */
#if NDIMS == 3
kh = kn / 2
kh = ni / 2
#endif /* NDIMS == 3 */
#ifdef MPI
@ -1331,8 +1331,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : in , jn , kn
use coordinates , only : ni => ncells, ng => nghosts
use coordinates , only : im , jm , km
use coordinates , only : faces_gc, faces_dc
use equations , only : nv
@ -1385,9 +1384,9 @@ module boundaries
! calculate half sizes
!
ih = in / 2
jh = jn / 2
kh = kn / 2
ih = ni / 2
jh = ni / 2
kh = ni / 2
#ifdef MPI
! prepare the block exchange structures
@ -1708,8 +1707,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : in , jn , kn
use coordinates , only : ni => ncells, ng => nghosts
use coordinates , only : im , jm , km
use coordinates , only : faces_gr
use equations , only : nv
@ -1760,9 +1758,9 @@ module boundaries
! calculate half sizes
!
ih = in / 2
jh = jn / 2
kh = kn / 2
ih = ni / 2
jh = ni / 2
kh = ni / 2
#ifdef MPI
! prepare the block exchange structures
@ -2074,8 +2072,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : in , jn , kn
use coordinates , only : ni => ncells, ng => nghosts
use coordinates , only : im , jm , km
use coordinates , only : faces_gp
use equations , only : nv
@ -2127,9 +2124,9 @@ module boundaries
! calculate the sizes
!
ih = in + ng
jh = jn + ng
kh = kn + ng
ih = ni + ng
jh = ni + ng
kh = ni + ng
#ifdef MPI
! prepare the block exchange structures
@ -2505,8 +2502,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : in, jn, kn
use coordinates , only : ni => ncells, ng => nghosts
use coordinates , only : im, jm, km
use coordinates , only : edges_gc, edges_dc
use equations , only : nv
@ -2559,10 +2555,10 @@ module boundaries
! calculate half sizes
!
ih = in / 2
jh = jn / 2
ih = ni / 2
jh = ni / 2
#if NDIMS == 3
kh = kn / 2
kh = ni / 2
#endif /* NDIMS == 3 */
#ifdef MPI
@ -2967,8 +2963,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : in , jn , kn
use coordinates , only : ni => ncells, ng => nghosts
use coordinates , only : im , jm , km
use coordinates , only : edges_gr
use equations , only : nv
@ -3019,10 +3014,10 @@ module boundaries
! calculate half sizes
!
ih = in / 2
jh = jn / 2
ih = ni / 2
jh = ni / 2
#if NDIMS == 3
kh = kn / 2
kh = ni / 2
#endif /* NDIMS == 3 */
#ifdef MPI
@ -3410,8 +3405,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : in , jn , kn
use coordinates , only : ni => ncells, ng => nghosts
use coordinates , only : im , jm , km
use coordinates , only : edges_gp
use equations , only : nv
@ -3463,10 +3457,10 @@ module boundaries
! calculate the sizes
!
ih = in + ng
jh = jn + ng
ih = ni + ng
jh = ni + ng
#if NDIMS == 3
kh = kn + ng
kh = ni + ng
#endif /* NDIMS == 3 */
#ifdef MPI
@ -3925,7 +3919,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : ng => nghosts
use coordinates , only : im, jm, km
use coordinates , only : corners_gc, corners_dc
use equations , only : nv
@ -4323,7 +4317,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : ng => nghosts
use coordinates , only : im , jm , km
use coordinates , only : corners_gr
use equations , only : nv
@ -4696,7 +4690,7 @@ module boundaries
use blocks , only : block_meta, block_data, block_leaf
use blocks , only : list_meta, list_leaf
use blocks , only : block_info, pointer_info
use coordinates , only : ng
use coordinates , only : ng => nghosts
use coordinates , only : im , jm , km
use coordinates , only : corners_gp
use equations , only : nv
@ -5078,7 +5072,7 @@ module boundaries
! import external procedures and variables
!
use coordinates , only : im , jm , km , ng
use coordinates , only : im , jm , km , ng => nghosts
use coordinates , only : ib , jb , kb , ie , je , ke
use coordinates , only : ibl, jbl, kbl, ieu, jeu, keu
use equations , only : nv
@ -5691,9 +5685,7 @@ module boundaries
! import external procedures and variables
!
use coordinates , only : ng, nd
use coordinates , only : in
use coordinates , only : in , jn , kn
use coordinates , only : ni => ncells, ng => nghosts, nd
use coordinates , only : im , jm , km
use coordinates , only : faces_dr
use equations , only : nv
@ -5736,8 +5728,8 @@ module boundaries
! calculate half sizes
!
jh = jn / 2
kh = kn / 2
jh = ni / 2
kh = ni / 2
! restrict the face region to the output array
!
@ -5755,8 +5747,8 @@ module boundaries
! calculate half sizes
!
ih = in / 2
kh = kn / 2
ih = ni / 2
kh = ni / 2
! restrict the face region to the output array
!
@ -5774,8 +5766,8 @@ module boundaries
! calculate half sizes
!
ih = in / 2
jh = jn / 2
ih = ni / 2
jh = ni / 2
! restrict the face region to the output array
!
@ -5816,9 +5808,7 @@ module boundaries
! import external procedures and variables
!
use coordinates , only : ng, nh
use coordinates , only : in
use coordinates , only : in , jn , kn
use coordinates , only : nh
use coordinates , only : im , jm , km
use coordinates , only : faces_dp
use equations , only : nv, idn, ipr
@ -5973,9 +5963,7 @@ module boundaries
! import external procedures and variables
!
use coordinates , only : ng, nd
use coordinates , only : in
use coordinates , only : in , jn , kn
use coordinates , only : ni => ncells, ng => nghosts, nd
use coordinates , only : im , jm , km
use coordinates , only : edges_dr
use equations , only : nv
@ -6028,7 +6016,7 @@ module boundaries
! calculate half size
!
ih = in / 2
ih = ni / 2
! restrict the edge region to the output array
!
@ -6055,7 +6043,7 @@ module boundaries
! calculate half size
!
jh = jn / 2
jh = ni / 2
! restrict the edge region to the output array
!
@ -6083,7 +6071,7 @@ module boundaries
! calculate half size
!
kh = kn / 2
kh = ni / 2
! restrict the edge region to the output array
!
@ -6125,9 +6113,7 @@ module boundaries
! import external procedures and variables
!
use coordinates , only : ng, nh
use coordinates , only : in
use coordinates , only : in , jn , kn
use coordinates , only : nh
use coordinates , only : im , jm , km
use coordinates , only : edges_dp
use equations , only : nv, idn, ipr
@ -6310,7 +6296,7 @@ module boundaries
! import external procedures and variables
!
use coordinates , only : ng, nd
use coordinates , only : ng => nghosts, nd
use coordinates , only : im , jm , km
use coordinates , only : corners_dr
use equations , only : nv
@ -6405,7 +6391,7 @@ module boundaries
! import external procedures and variables
!
use coordinates , only : ng, nh
use coordinates , only : ng => nghosts, nh
use coordinates , only : im , jm , km
use coordinates , only : corners_dp
use equations , only : nv, idn, ipr
@ -6594,7 +6580,7 @@ module boundaries
! import external procedures and variables
!
use blocks , only : block_data
use coordinates , only : in, jn, kn
use coordinates , only : ni => ncells
use equations , only : nv
! local variables are not implicit by default
@ -6621,13 +6607,13 @@ module boundaries
#if NDIMS == 2
! average fluxes from higher level neighbor
!
fb(1:nv,:,:) = (fn(1:nv,1:jn:2,1:kn) + fn(1:nv,2:jn:2,1:kn)) / 2.0d+00
fb(1:nv,:,:) = (fn(1:nv,1:ni:2, : ) + fn(1:nv,2:ni:2, : )) / 2.0d+00
#endif /* NDIMS == 2 */
#if NDIMS == 3
! average fluxes from higher level neighbor
!
fb(1:nv,:,:) = ((fn(1:nv,1:in:2,1:kn:2) + fn(1:nv,2:in:2,2:kn:2)) &
+ (fn(1:nv,1:in:2,2:kn:2) + fn(1:nv,2:in:2,1:kn:2))) &
fb(1:nv,:,:) = ((fn(1:nv,1:ni:2,1:ni:2) + fn(1:nv,2:ni:2,2:ni:2)) &
+ (fn(1:nv,1:ni:2,2:ni:2) + fn(1:nv,2:ni:2,1:ni:2))) &
/ 4.0d+00
#endif /* NDIMS == 3 */
@ -6638,13 +6624,13 @@ module boundaries
#if NDIMS == 2
! average fluxes from higher level neighbor
!
fb(1:nv,:,:) = (fn(1:nv,1:in:2,1:kn) + fn(1:nv,2:in:2,1:kn)) / 2.0d+00
fb(1:nv,:,:) = (fn(1:nv,1:ni:2, : ) + fn(1:nv,2:ni:2, : )) / 2.0d+00
#endif /* NDIMS == 2 */
#if NDIMS == 3
! average fluxes from higher level neighbor
!
fb(1:nv,:,:) = ((fn(1:nv,1:in:2,1:kn:2) + fn(1:nv,2:in:2,2:kn:2)) &
+ (fn(1:nv,1:in:2,2:kn:2) + fn(1:nv,2:in:2,1:kn:2))) &
fb(1:nv,:,:) = ((fn(1:nv,1:ni:2,1:ni:2) + fn(1:nv,2:ni:2,2:ni:2)) &
+ (fn(1:nv,1:ni:2,2:ni:2) + fn(1:nv,2:ni:2,1:ni:2))) &
/ 4.0d+00
#endif /* NDIMS == 3 */
@ -6655,8 +6641,8 @@ module boundaries
! average fluxes from higher level neighbor
!
fb(1:nv,:,:) = ((fn(1:nv,1:in:2,1:jn:2) + fn(1:nv,2:in:2,2:jn:2)) &
+ (fn(1:nv,1:in:2,2:jn:2) + fn(1:nv,2:in:2,1:jn:2))) &
fb(1:nv,:,:) = ((fn(1:nv,1:ni:2,1:ni:2) + fn(1:nv,2:ni:2,2:ni:2)) &
+ (fn(1:nv,1:ni:2,2:ni:2) + fn(1:nv,2:ni:2,1:ni:2))) &
/ 4.0d+00
#endif /* NDIMS == 3 */
@ -6688,7 +6674,7 @@ module boundaries
! include external variables
!
use blocks , only : block_data, list_data
use coordinates , only : im , jm , km , in , jn , kn
use coordinates , only : ni => ncells, im , jm , km
use coordinates , only : ib , jb , kb , ie , je , ke
use coordinates , only : ibl, jbl, kbl, ieu, jeu, keu
use equations , only : nv
@ -6741,13 +6727,13 @@ module boundaries
! update remaining left layers of the X boundary
!
do i = 1, ibl
call prim2cons(jn, pdata%q(1:nv,i,jb:je,k), pdata%u(1:nv,i,jb:je,k))
call prim2cons(ni, pdata%q(1:nv,i,jb:je,k), pdata%u(1:nv,i,jb:je,k))
end do ! i = 1, ibl
! update remaining right layers of the X boundary
!
do i = ieu, im
call prim2cons(jn, pdata%q(1:nv,i,jb:je,k), pdata%u(1:nv,i,jb:je,k))
call prim2cons(ni, pdata%q(1:nv,i,jb:je,k), pdata%u(1:nv,i,jb:je,k))
end do ! i = 1, ibl
end do ! k = 1, km
@ -6760,13 +6746,13 @@ module boundaries
! update the remaining front layers of the Z boundary
!
do k = 1, kbl
call prim2cons(in, pdata%q(1:nv,ib:ie,j,k), pdata%u(1:nv,ib:ie,j,k))
call prim2cons(ni, pdata%q(1:nv,ib:ie,j,k), pdata%u(1:nv,ib:ie,j,k))
end do ! k = 1, kbl
! update the remaining back layers of the Z boundary
!
do k = keu, km
call prim2cons(in, pdata%q(1:nv,ib:ie,j,k), pdata%u(1:nv,ib:ie,j,k))
call prim2cons(ni, pdata%q(1:nv,ib:ie,j,k), pdata%u(1:nv,ib:ie,j,k))
end do ! k = keu, km
end do ! j = jb, je

View File

@ -37,13 +37,13 @@ module coordinates
! MODULE PARAMETERS:
! =================
!
! the domain block dimensions
! the number of interior cells along each block dimension
!
integer , save :: nc = 8, in = 8, jn = 8, kn = 1
integer , save :: ncells = 8
! the number of ghost zones
! the number of ghost zones at each side of the block
!
integer , save :: ng = 2, nh = 1, nd = 4
integer , save :: nghosts = 2, nh = 1, nd = 4
! the domain block dimensions including the ghost zones
!
@ -152,8 +152,8 @@ module coordinates
!
! Arguments:
!
! ncells - the number of cells per dimension in a block;
! nghosts - the number of ghosts cells at each side in a block;
! ni - the number of interior cells per dimension in a block;
! ng - the number of ghosts cells at each side in a block;
! tlevel - the maximum refinement level for coordinate generation;
! bdims - the base level dimension in blocks;
! xmn, xmx, ymn, ymx, zmn, zmx -
@ -163,7 +163,7 @@ module coordinates
!
!===============================================================================
!
subroutine initialize_coordinates(ncells, nghosts, tlevel, bdims, &
subroutine initialize_coordinates(ni, ng, tlevel, bdims, &
xmn, xmx, ymn, ymx, zmn, zmx, verbose, iret)
! include external procedures
@ -177,7 +177,7 @@ module coordinates
! subroutine arguments
!
logical , intent(in) :: verbose
integer , intent(in) :: ncells, nghosts, tlevel
integer , intent(in) :: ni, ng, tlevel
integer, dimension(3), intent(in) :: bdims
real(kind=8) , intent(in) :: xmn, xmx, ymn, ymx, zmn, zmx
integer , intent(inout) :: iret
@ -185,9 +185,9 @@ module coordinates
! local variables
!
character(len=32) :: lbndry, ubndry
integer :: i, j, k, l, p, q, r, ff
integer :: i, j, k, l, p, q, r
integer :: fi, fj, fk
integer :: ni, nj, nk, nm, np, nr, nq, ns, nt, nu
integer :: nn, nm, np, nr, nq, ns, nt, nu
logical :: info
!
!-------------------------------------------------------------------------------
@ -195,16 +195,8 @@ module coordinates
! set the number of cells along each block dimension and
! the number of ghost zones
!
nc = ncells
ng = nghosts
! set the block dimensions
!
in = nc
jn = nc
#if NDIMS == 3
kn = nc
#endif /* NDIMS == 3 */
ncells = ni
nghosts = ng
! calculate half and double of the number of ghose zones
!
@ -213,10 +205,10 @@ module coordinates
! calculate the block dimensions including ghost cells
!
im = in + 2 * ng
jm = jn + 2 * ng
im = ni + 2 * ng
jm = ni + 2 * ng
#if NDIMS == 3
km = kn + 2 * ng
km = ni + 2 * ng
#endif /* NDIMS == 3 */
! prepare indices
@ -231,14 +223,14 @@ module coordinates
!
ih = im / 2
ib = ng + 1
ie = ng + in
ie = ng + ni
ibl = ib - 1
ibu = ib + ng - 1
iel = ie - ng + 1
ieu = ie + 1
jh = jm / 2
jb = ng + 1
je = ng + jn
je = ng + ni
jbl = jb - 1
jbu = jb + ng - 1
jel = je - ng + 1
@ -246,7 +238,7 @@ module coordinates
#if NDIMS == 3
kh = km / 2
kb = ng + 1
ke = ng + kn
ke = ng + ni
kbl = kb - 1
kbu = kb + ng - 1
kel = ke - ng + 1
@ -350,17 +342,14 @@ module coordinates
! calculate the block resolution at each level
!
ff = 2**(l - 1)
ni = in * ff
nj = jn * ff
nk = kn * ff
nn = ni * 2**(l - 1)
! calculate the cell sizes for each level
!
adx (l) = xlen / (domain_base_dims(1) * ni)
ady (l) = ylen / (domain_base_dims(2) * nj)
adx (l) = xlen / (domain_base_dims(1) * nn)
ady (l) = ylen / (domain_base_dims(2) * nn)
#if NDIMS == 3
adz (l) = zlen / (domain_base_dims(3) * nk)
adz (l) = zlen / (domain_base_dims(3) * nn)
#endif /* NDIMS == 3 */
#if NDIMS == 2
adr (l) = sqrt(adx(l)**2 + ady(l)**2)
@ -393,11 +382,11 @@ module coordinates
! initialize ghost subarray indices
!
np = nc + ng
nm = nc - ng
nr = nc - nd
nq = nc - nh
ns = nc / 2
np = ni + ng
nm = ni - ng
nr = ni - nd
nq = ni - nh
ns = ni / 2
nt = np / 2
nu = nm / 2
#if NDIMS == 2
@ -441,8 +430,8 @@ module coordinates
edges_gr(i,j,1)%u(:) = edges_gr(i,j,1)%l(:) + (/ ns, ng /) - 1
edges_gr(i,j,2)%u(:) = edges_gr(i,j,2)%l(:) + (/ ng, ns /) - 1
edges_dr(i,q,1)%u(:) = edges_dr(i,q,1)%l(:) + (/ nc, nd /) - 1
edges_dr(p,j,2)%u(:) = edges_dr(p,j,2)%l(:) + (/ nd, nc /) - 1
edges_dr(i,q,1)%u(:) = edges_dr(i,q,1)%l(:) + (/ ni, nd /) - 1
edges_dr(p,j,2)%u(:) = edges_dr(p,j,2)%l(:) + (/ nd, ni /) - 1
! for edge prolongation
!
@ -568,9 +557,9 @@ module coordinates
faces_gr(i,j,k,2)%u(:) = faces_gr(i,j,k,2)%l(:) + (/ ns, ng, ns /) - 1
faces_gr(i,j,k,3)%u(:) = faces_gr(i,j,k,3)%l(:) + (/ ns, ns, ng /) - 1
faces_dr(p,j,k,1)%u(:) = faces_dr(p,j,k,1)%l(:) + (/ nd, nc, nc /) - 1
faces_dr(i,q,k,2)%u(:) = faces_dr(i,q,k,2)%l(:) + (/ nc, nd, nc /) - 1
faces_dr(i,j,r,3)%u(:) = faces_dr(i,j,r,3)%l(:) + (/ nc, nc, nd /) - 1
faces_dr(p,j,k,1)%u(:) = faces_dr(p,j,k,1)%l(:) + (/ nd, ni, ni /) - 1
faces_dr(i,q,k,2)%u(:) = faces_dr(i,q,k,2)%l(:) + (/ ni, nd, ni /) - 1
faces_dr(i,j,r,3)%u(:) = faces_dr(i,j,r,3)%l(:) + (/ ni, ni, nd /) - 1
! for face prolongation
!
@ -658,9 +647,9 @@ module coordinates
edges_gr(i,j,k,2)%u(:) = edges_gr(i,j,k,2)%l(:) + (/ ng, ns, ng /) - 1
edges_gr(i,j,k,3)%u(:) = edges_gr(i,j,k,3)%l(:) + (/ ng, ng, ns /) - 1
edges_dr(i,q,r,1)%u(:) = edges_dr(i,q,r,1)%l(:) + (/ nc, nd, nd /) - 1
edges_dr(p,j,r,2)%u(:) = edges_dr(p,j,r,2)%l(:) + (/ nd, nc, nd /) - 1
edges_dr(p,q,k,3)%u(:) = edges_dr(p,q,k,3)%l(:) + (/ nd, nd, nc /) - 1
edges_dr(i,q,r,1)%u(:) = edges_dr(i,q,r,1)%l(:) + (/ ni, nd, nd /) - 1
edges_dr(p,j,r,2)%u(:) = edges_dr(p,j,r,2)%l(:) + (/ nd, ni, nd /) - 1
edges_dr(p,q,k,3)%u(:) = edges_dr(p,q,k,3)%l(:) + (/ nd, nd, ni /) - 1
! for edge prolongation
!
@ -831,13 +820,13 @@ module coordinates
!
bm(:) = domain_base_dims(:)
rm(:) = bm(:) * 2**(maxlev - 1)
cm(:) = bm(:) * nc
fm(:) = rm(:) * nc
cm(:) = bm(:) * ncells
fm(:) = rm(:) * ncells
call print_section(verbose, "Geometry")
call print_parameter(verbose, "refinement to level" , maxlev)
call print_parameter(verbose, "number of block cells", nc )
call print_parameter(verbose, "number of ghost zones", ng )
call print_parameter(verbose, "refinement to level" , maxlev )
call print_parameter(verbose, "number of block cells", ncells )
call print_parameter(verbose, "number of ghost zones", nghosts)
write(msg,"(i0)") maxval(fm(1:NDIMS))
p = len(trim(adjustl(msg)))
write(msg,"(i0)") maxval(rm(1:NDIMS))

View File

@ -1150,7 +1150,7 @@ module equations
! include external procedures and variables
!
use coordinates, only : im , jm , km , in , jn , kn
use coordinates, only : ncells, im , jm , km
use coordinates, only : ib , jb , kb , ie , je , ke
use coordinates, only : ibl, jbl, kbl, ieu, jeu, keu
@ -1176,7 +1176,7 @@ module equations
! convert conserved variables to primitive ones
!
call cons2prim(in, uu(1:nv,ib:ie,j,k), qq(1:nv,ib:ie,j,k))
call cons2prim(ncells, uu(1:nv,ib:ie,j,k), qq(1:nv,ib:ie,j,k))
end do ! j = jb, je
end do ! k = kb, ke

View File

@ -256,7 +256,7 @@ module integrals
! import external variables and subroutines
!
use blocks , only : block_meta, block_data, list_data
use coordinates , only : in, jn, kn, ib, jb, kb, ie, je, ke
use coordinates , only : ni => ncells, ib, jb, kb, ie, je, ke
use coordinates , only : advol, voli
use equations , only : idn, ipr, ivx, ivy, ivz, ibx, iby, ibz, ibp
use equations , only : ien, imx, imy, imz
@ -289,7 +289,11 @@ module integrals
! local arrays
!
real(kind=8), dimension(narr) :: inarr, avarr, mnarr, mxarr
real(kind=8), dimension(in,jn,kn) :: vel, mag, sqd, tmp
#if NDIMS == 3
real(kind=8), dimension(ni,ni,ni) :: vel, mag, sqd, tmp
#else /* NDIMS == 3 */
real(kind=8), dimension(ni,ni, 1) :: vel, mag, sqd, tmp
#endif /* NDIMS == 3 */
! parameters
!

View File

@ -1773,7 +1773,7 @@ module io
use blocks , only : get_mblocks, get_dblocks, get_nleafs
use blocks , only : get_last_id
use coordinates , only : minlev, maxlev
use coordinates , only : nc, ng, in, jn, kn
use coordinates , only : ncells, nghosts
use coordinates , only : ddims => domain_base_dims
use coordinates , only : xmin, xmax, ymin, ymax, zmin, zmax
use coordinates , only : periodic
@ -1801,6 +1801,7 @@ module io
! local vectors
!
integer, dimension(3) :: bdims = 1
integer, dimension(3) :: per
! local allocatable arrays
@ -1854,8 +1855,8 @@ module io
call write_attribute(gid, 'mblocks' , get_mblocks())
call write_attribute(gid, 'dblocks' , get_dblocks())
call write_attribute(gid, 'nleafs' , get_nleafs() )
call write_attribute(gid, 'ncells' , nc )
call write_attribute(gid, 'nghosts' , ng )
call write_attribute(gid, 'ncells' , ncells )
call write_attribute(gid, 'nghosts' , nghosts )
call write_attribute(gid, 'minlev' , minlev )
call write_attribute(gid, 'maxlev' , maxlev )
call write_attribute(gid, 'nprocs' , nprocs )
@ -1885,7 +1886,8 @@ module io
! store the vector attributes
!
call write_attribute(gid, 'dims' , (/ in, jn, kn /))
bdims(1:NDIMS) = ncells
call write_attribute(gid, 'dims' , bdims)
call write_attribute(gid, 'domain_base_dims', ddims)
! store random number generator seed values
@ -1951,7 +1953,7 @@ module io
use blocks , only : append_metablock
use blocks , only : set_last_id, get_last_id
use blocks , only : get_mblocks, get_dblocks, get_nleafs
use coordinates , only : nc, ng, in, jn, kn
use coordinates , only : ncells, nghosts
use coordinates , only : xmin, xmax, ymin, ymax, zmin, zmax
use evolution , only : step, time, dt, dtn
use hdf5 , only : hid_t, hsize_t
@ -2037,14 +2039,14 @@ module io
! check the block dimensions
!
if (lncells /= nc) then
if (lncells /= ncells) then
write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "The block dimensions do not match!"
end if
! check the number of ghost layers
!
if (lnghost /= ng) then
if (lnghost /= nghosts) then
write(error_unit,"('[',a,']: ',a)") trim(loc) &
, "The number of ghost layers does not match!"
end if
@ -3303,7 +3305,7 @@ module io
!
use blocks , only : block_data, list_data
use blocks , only : get_dblocks
use coordinates , only : im, jm, km, in, jn, kn
use coordinates , only : ni => ncells, im, jm, km
use coordinates , only : ib, jb, kb, ie, je, ke
use equations , only : nv, pvars
use hdf5 , only : hid_t, hsize_t
@ -3321,7 +3323,7 @@ module io
! HDF5 variables
!
integer(hid_t) :: gid
integer(hsize_t) :: dm(4)
integer(hsize_t) :: dm(4) = 1
! local variables
!
@ -3371,9 +3373,11 @@ module io
ju = jm
ku = km
else
dm(2) = in
dm(3) = jn
dm(4) = kn
dm(2) = ni
dm(3) = ni
#if NDIMS == 3
dm(4) = ni
#endif /* NDIMS == 3 */
il = ib
jl = jb
@ -3479,7 +3483,7 @@ module io
!
use blocks , only : block_data, list_data
use blocks , only : get_dblocks
use coordinates , only : im, jm, km, in, jn, kn
use coordinates , only : ni => ncells, im, jm, km
use coordinates , only : ib, jb, kb, ie, je, ke
use equations , only : nv, cvars
use hdf5 , only : hid_t, hsize_t
@ -3497,7 +3501,7 @@ module io
! HDF5 variables
!
integer(hid_t) :: gid
integer(hsize_t) :: dm(4)
integer(hsize_t) :: dm(4) = 1
! local variables
!
@ -3547,9 +3551,11 @@ module io
ju = jm
ku = km
else
dm(2) = in
dm(3) = jn
dm(4) = kn
dm(2) = ni
dm(3) = ni
#if NDIMS == 3
dm(4) = ni
#endif /* NDIMS == 3 */
il = ib
jl = jb
@ -7174,9 +7180,8 @@ module io
use blocks , only : get_dblocks
use equations , only : nv, pvars
use mpitools , only : nproc
use coordinates , only : in, jn, kn
use coordinates , only : ni => ncells, ng => nghosts
use coordinates , only : adx, ady, adz
use coordinates , only : ng
use evolution , only : time
! local variables are not implicit by default
@ -7232,27 +7237,27 @@ module io
! prepare dimensions
!
ip = in + 1
jp = jn + 1
ip = ni + 1
jp = ni + 1
#if NDIMS == 3
kp = kn + 1
kp = ni + 1
#endif /* NDIMS == 3 */
#if NDIMS == 3
write(stmp, "(1i8)") kn
write(ttmp, "(1i8)") jn
write(stmp, "(1i8)") ni
write(ttmp, "(1i8)") ni
stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
write(ttmp, "(1i8)") in
write(ttmp, "(1i8)") ni
#else /* NDIMS == 3 */
write(stmp, "(1i8)") jn
write(ttmp, "(1i8)") in
write(stmp, "(1i8)") ni
write(ttmp, "(1i8)") ni
#endif /* NDIMS == 3 */
bdim = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
write(stmp, "(1i8)") kn
write(ttmp, "(1i8)") jn
write(stmp, "(1i8)") ni
write(ttmp, "(1i8)") ni
stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
write(ttmp, "(1i8)") in
write(ttmp, "(1i8)") ni
stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
write(ttmp, "(1i8)") get_dblocks()
sdim = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
@ -7261,15 +7266,15 @@ module io
!
#if NDIMS == 3
if (with_ghosts) then
slab(:) = (/ ng, ng, ng, -1, 1, 1, 1, 1, kn, jn, in, 1 /)
slab(:) = (/ ng, ng, ng, -1, 1, 1, 1, 1, ni, ni, ni, 1 /)
else
slab(:) = (/ 0, 0, 0, -1, 1, 1, 1, 1, kn, jn, in, 1 /)
slab(:) = (/ 0, 0, 0, -1, 1, 1, 1, 1, ni, ni, ni, 1 /)
end if
#else /* NDIMS == 3 */
if (with_ghosts) then
slab(:) = (/ 0, ng, ng, -1, 1, 1, 1, 1, 1, jn, in, 1 /)
slab(:) = (/ 0, ng, ng, -1, 1, 1, 1, 1, 1, ni, ni, 1 /)
else
slab(:) = (/ 0, 0, 0, -1, 1, 1, 1, 1, 1, jn, in, 1 /)
slab(:) = (/ 0, 0, 0, -1, 1, 1, 1, 1, 1, ni, ni, 1 /)
end if
#endif /* NDIMS == 3 */

View File

@ -261,7 +261,7 @@ module mesh
use blocks , only : ndims
use blocks , only : get_mblocks, get_nleafs
use coordinates , only : ddims => domain_base_dims
use coordinates , only : ng, nd, in, jn, kn, im, jm, km
use coordinates , only : ni => ncells, ng => nghosts, nd, im, jm, km
use coordinates , only : toplev
use mpitools , only : master, nprocs
@ -338,10 +338,10 @@ module mesh
!
ff = 2**(toplev - 1)
fcv = 1.0d+00 / n
fef = 1.0d+00 * (ddims(1) * in * ff + nd) / im
fef = fef * (ddims(2) * jn * ff + nd) / jm
fef = 1.0d+00 * (ddims(1) * ni * ff + nd) / im
fef = fef * (ddims(2) * ni * ff + nd) / jm
#if NDIMS == 3
fef = fef * (ddims(3) * kn * ff + nd) / km
fef = fef * (ddims(3) * ni * ff + nd) / km
#endif /* NDIMS == 3 */
! reset the first execution flag
@ -1034,7 +1034,7 @@ module mesh
! import external procedures and variables
!
use blocks , only : block_meta, block_data, nchildren
use coordinates , only : ng, nh, in, jn, kn, im, jm, km
use coordinates , only : ncells, nghosts, nh, im, jm, km
use coordinates , only : ib, ie, jb, je, kb, ke
use equations , only : nv, idn, ien
use interpolations , only : limiter_prol
@ -1062,7 +1062,7 @@ module mesh
! local arrays
!
integer , dimension(3) :: dm
integer , dimension(3) :: dm = 1
real(kind=8), dimension(NDIMS) :: du
! local allocatable arrays
@ -1087,10 +1087,7 @@ module mesh
! prepare dimensions
!
dm(:) = 2 * ((/ in, jn, kn /) + ng)
#if NDIMS == 2
dm(3) = 1
#endif /* NDIMS == 2 */
dm(1:NDIMS) = 2 * (ncells + nghosts)
! allocate array for the expanded arrays
!
@ -1203,10 +1200,10 @@ module mesh
! calculate indices of the current child subdomain
!
il = 1 + ic * in
jl = 1 + jc * jn
il = 1 + ic * ncells
jl = 1 + jc * ncells
#if NDIMS == 3
kl = 1 + kc * kn
kl = 1 + kc * ncells
#endif /* NDIMS == 3 */
iu = il + im - 1
@ -1260,9 +1257,8 @@ module mesh
! import external procedures and variables
!
use blocks , only : ndims
use blocks , only : block_meta, block_data, nchildren
use coordinates , only : ng, nh, in, jn, kn, im, jm, km
use coordinates , only : nh, im, jm, km
use coordinates , only : ih, jh, kh, ib, jb, kb, ie, je, ke
use equations , only : nv
@ -1282,7 +1278,7 @@ module mesh
! local arrays
!
integer, dimension(ndims) :: pos
integer, dimension(NDIMS) :: pos
! local pointers
!
@ -1310,7 +1306,7 @@ module mesh
! obtain the child position in the parent block
!
pos(1:ndims) = pchild%meta%pos(1:ndims)
pos(1:NDIMS) = pchild%meta%pos(1:NDIMS)
! calculate the bound indices of the source nad destination arrays
!