MESH: Move mesh statistics to be handled in INTEGRALS.
Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
This commit is contained in:
parent
1979eae289
commit
dd3da22460
@ -41,7 +41,7 @@ program amun
|
|||||||
use io , only : read_restart_snapshot, write_restart_snapshot
|
use io , only : read_restart_snapshot, write_restart_snapshot
|
||||||
use io , only : write_snapshot, update_dtp
|
use io , only : write_snapshot, update_dtp
|
||||||
use iso_fortran_env, only : error_unit
|
use iso_fortran_env, only : error_unit
|
||||||
use mesh , only : generate_mesh, store_mesh_stats
|
use mesh , only : generate_mesh
|
||||||
use mpitools , only : initialize_mpitools, finalize_mpitools
|
use mpitools , only : initialize_mpitools, finalize_mpitools
|
||||||
#ifdef MPI
|
#ifdef MPI
|
||||||
use mpitools , only : reduce_sum
|
use mpitools , only : reduce_sum
|
||||||
@ -237,8 +237,6 @@ program amun
|
|||||||
call user_time_statistics(time)
|
call user_time_statistics(time)
|
||||||
call write_snapshot(name)
|
call write_snapshot(name)
|
||||||
|
|
||||||
call store_mesh_stats(step, time)
|
|
||||||
|
|
||||||
end if
|
end if
|
||||||
|
|
||||||
call stop_timer(iin)
|
call stop_timer(iin)
|
||||||
@ -305,10 +303,6 @@ program amun
|
|||||||
!
|
!
|
||||||
call new_time_step()
|
call new_time_step()
|
||||||
|
|
||||||
! store the mesh statistics
|
|
||||||
!
|
|
||||||
call store_mesh_stats(step, time)
|
|
||||||
|
|
||||||
! store the time statistics and the run snapshot
|
! store the time statistics and the run snapshot
|
||||||
!
|
!
|
||||||
call store_integrals()
|
call store_integrals()
|
||||||
|
@ -35,30 +35,31 @@ module integrals
|
|||||||
! MODULE PARAMETERS:
|
! MODULE PARAMETERS:
|
||||||
! =================
|
! =================
|
||||||
!
|
!
|
||||||
! funit - a file handler to the integrals file;
|
! munit - a file handler for the mesh statistics file;
|
||||||
! sunit - a file handler to the statistics file;
|
! funit - a file handler for the integrals file;
|
||||||
! eunit - a file handler to the errors file;
|
! sunit - a file handler for the statistics file;
|
||||||
! iintd - the number of steps between subsequent intervals storing;
|
! eunit - a file handler for the errors file;
|
||||||
!
|
!
|
||||||
|
integer(kind=4), save :: munit = 10
|
||||||
integer(kind=4), save :: funit = 11
|
integer(kind=4), save :: funit = 11
|
||||||
integer(kind=4), save :: sunit = 12
|
integer(kind=4), save :: sunit = 12
|
||||||
integer(kind=4), save :: eunit = 13
|
integer(kind=4), save :: eunit = 13
|
||||||
|
|
||||||
|
! iintd - the number of steps between subsequent intervals;
|
||||||
|
!
|
||||||
integer(kind=4), save :: iintd = 1
|
integer(kind=4), save :: iintd = 1
|
||||||
|
|
||||||
! flag indicating if the files are actually written to disk
|
! the mesh coverage and efficiency factors
|
||||||
!
|
!
|
||||||
logical, save :: stored = .false.
|
real(kind=8), save :: fcv = 1.0d+00, fef = 1.0d+00
|
||||||
|
|
||||||
! the format string
|
! the format strings
|
||||||
!
|
!
|
||||||
|
character(len=32), save :: mfmt1, mfmt2, mfmt3
|
||||||
character(len=32), save :: efmt
|
character(len=32), save :: efmt
|
||||||
|
|
||||||
! by default everything is private
|
|
||||||
!
|
|
||||||
private
|
private
|
||||||
|
|
||||||
! declare public subroutines
|
|
||||||
!
|
|
||||||
public :: initialize_integrals, finalize_integrals
|
public :: initialize_integrals, finalize_integrals
|
||||||
public :: store_integrals
|
public :: store_integrals
|
||||||
|
|
||||||
@ -81,54 +82,91 @@ module integrals
|
|||||||
!
|
!
|
||||||
! Arguments:
|
! Arguments:
|
||||||
!
|
!
|
||||||
! store - flag indicating if the files should be actually written to disk;
|
! verbose - the verbose flag;
|
||||||
! irun - job execution counter;
|
! nrun - the number of resumed run;
|
||||||
! status - return flag of the procedure execution status;
|
! status - return flag of the procedure execution status;
|
||||||
!
|
!
|
||||||
!===============================================================================
|
!===============================================================================
|
||||||
!
|
!
|
||||||
subroutine initialize_integrals(store, irun, status)
|
subroutine initialize_integrals(verbose, nrun, status)
|
||||||
|
|
||||||
! import external variables and subroutines
|
use coordinates, only : toplev, ncells, bcells, nghosts, domain_base_dims
|
||||||
!
|
use equations , only : pvars, nf
|
||||||
use equations , only : pvars, nf
|
use evolution , only : error_control
|
||||||
use evolution , only : error_control
|
use mpitools , only : master, nprocs
|
||||||
use parameters, only : get_parameter
|
use parameters , only : get_parameter
|
||||||
|
|
||||||
! local variables are not implicit by default
|
|
||||||
!
|
|
||||||
implicit none
|
implicit none
|
||||||
|
|
||||||
! subroutine arguments
|
logical, intent(in) :: verbose
|
||||||
!
|
integer, intent(in) :: nrun
|
||||||
logical, intent(in) :: store
|
integer, intent(out) :: status
|
||||||
integer, intent(inout) :: irun
|
|
||||||
integer, intent(out) :: status
|
|
||||||
|
|
||||||
! local variables
|
|
||||||
!
|
|
||||||
character(len=32) :: fname, append, stmp
|
character(len=32) :: fname, append, stmp
|
||||||
logical :: flag
|
logical :: flag
|
||||||
!
|
integer :: l
|
||||||
|
|
||||||
!-------------------------------------------------------------------------------
|
!-------------------------------------------------------------------------------
|
||||||
!
|
!
|
||||||
status = 0
|
status = 0
|
||||||
|
|
||||||
! set stored flag
|
|
||||||
!
|
|
||||||
stored = store
|
|
||||||
|
|
||||||
! get the integrals storing interval
|
|
||||||
!
|
|
||||||
call get_parameter("integrals_interval", iintd)
|
call get_parameter("integrals_interval", iintd)
|
||||||
|
|
||||||
! make sure storage interval is larger than zero
|
|
||||||
!
|
|
||||||
iintd = max(1, iintd)
|
iintd = max(1, iintd)
|
||||||
|
|
||||||
! handle the integral file if store is set
|
! only master process writes the files to the disk
|
||||||
!
|
!
|
||||||
if (stored) then
|
if (master) then
|
||||||
|
|
||||||
|
! prepare the mesh statistics file
|
||||||
|
!
|
||||||
|
write(fname,"('mesh-statistics_',i2.2,'.dat')") nrun
|
||||||
|
|
||||||
|
#ifdef __INTEL_COMPILER
|
||||||
|
open(newunit=munit, file=fname, form='formatted', &
|
||||||
|
status='replace', buffered='yes')
|
||||||
|
#else /* __INTEL_COMPILER */
|
||||||
|
open(newunit=munit, file=fname, form='formatted', status='replace')
|
||||||
|
#endif /* __INTEL_COMPILER */
|
||||||
|
|
||||||
|
! write the header of the mesh statistics file
|
||||||
|
!
|
||||||
|
write(munit,"('#',5x,'step',4x,'time',11x,'coverage',7x,'efficiency',"// &
|
||||||
|
"5x,'blocks',5x,'leafs',1x,'|',1x," // &
|
||||||
|
"'block distribution per level')", advance="no")
|
||||||
|
#ifdef MPI
|
||||||
|
write(stmp,"(a,i0,a)") "(", max(1, 10 * toplev - 28), "x,'|',1x,a)"
|
||||||
|
write(munit,stmp) "block distribution per process"
|
||||||
|
#endif /* MPI */
|
||||||
|
write(munit,"('#',75x,'|')", advance="no")
|
||||||
|
do l = 1, toplev
|
||||||
|
write(munit,"(1x,i9)", advance="no") l
|
||||||
|
end do
|
||||||
|
#ifdef MPI
|
||||||
|
write(stmp,"(a,i0,a)") "(", max(1, 10 * (3 - toplev)), "x,'|')"
|
||||||
|
write(munit,stmp, advance="no")
|
||||||
|
do l = 0, nprocs - 1
|
||||||
|
write(munit,"(1x,i9)", advance="no") l
|
||||||
|
end do
|
||||||
|
#endif /* MPI */
|
||||||
|
write(munit,"('')")
|
||||||
|
|
||||||
|
! calculate the coverage and efficiency factors
|
||||||
|
!
|
||||||
|
fcv = 1.0d+00 / (product(domain_base_dims) * 2**(NDIMS*(toplev - 1)))
|
||||||
|
fef = 1.0d+00
|
||||||
|
do l = 1, NDIMS
|
||||||
|
fef = fef * (domain_base_dims(l) * ncells * 2**(toplev - 1) &
|
||||||
|
+ 2 * nghosts) / bcells
|
||||||
|
end do
|
||||||
|
|
||||||
|
! prepare the format strings for mesh statistics
|
||||||
|
!
|
||||||
|
mfmt1 = "(1x,i9,3(1x,1es14.6e3),2(1x,i9))"
|
||||||
|
write (mfmt2,"(a,i0,a)") '(2x,', toplev, '(1x,i9))'
|
||||||
|
#ifdef MPI
|
||||||
|
write(stmp,"(a,i0,a)") "(", max(2, 10 * (3 - toplev) + 1), "x,"
|
||||||
|
write(mfmt3,"(a,i0,a)") trim(stmp), nprocs, '(1x,i9))'
|
||||||
|
#endif /* MPI */
|
||||||
|
|
||||||
! depending on the append parameter, choose the right file initialization for
|
! depending on the append parameter, choose the right file initialization for
|
||||||
! the integrals file
|
! the integrals file
|
||||||
@ -140,13 +178,13 @@ module integrals
|
|||||||
write(fname, "('integrals.dat')")
|
write(fname, "('integrals.dat')")
|
||||||
inquire(file=fname, exist=flag)
|
inquire(file=fname, exist=flag)
|
||||||
case default
|
case default
|
||||||
write(fname, "('integrals_',i2.2,'.dat')") irun
|
write(fname, "('integrals_',i2.2,'.dat')") nrun
|
||||||
flag = .false.
|
flag = .false.
|
||||||
end select
|
end select
|
||||||
|
|
||||||
! check if the file exists; if not, create a new one, otherwise move to the end
|
! check if the file exists; if not, create a new one, otherwise move to the end
|
||||||
!
|
!
|
||||||
if (flag .and. irun > 1) then
|
if (flag .and. nrun > 1) then
|
||||||
#ifdef __INTEL_COMPILER
|
#ifdef __INTEL_COMPILER
|
||||||
open(newunit=funit, file=fname, form='formatted', status='old', &
|
open(newunit=funit, file=fname, form='formatted', status='old', &
|
||||||
position='append', buffered='yes')
|
position='append', buffered='yes')
|
||||||
@ -182,13 +220,13 @@ module integrals
|
|||||||
write(fname, "('statistics.dat')")
|
write(fname, "('statistics.dat')")
|
||||||
inquire(file=fname, exist=flag)
|
inquire(file=fname, exist=flag)
|
||||||
case default
|
case default
|
||||||
write(fname, "('statistics_',i2.2,'.dat')") irun
|
write(fname, "('statistics_',i2.2,'.dat')") nrun
|
||||||
flag = .false.
|
flag = .false.
|
||||||
end select
|
end select
|
||||||
|
|
||||||
! check if the file exists; if not, create a new one, otherwise move to the end
|
! check if the file exists; if not, create a new one, otherwise move to the end
|
||||||
!
|
!
|
||||||
if (flag .and. irun > 1) then
|
if (flag .and. nrun > 1) then
|
||||||
#ifdef __INTEL_COMPILER
|
#ifdef __INTEL_COMPILER
|
||||||
open(newunit=sunit, file=fname, form='formatted', status='old', &
|
open(newunit=sunit, file=fname, form='formatted', status='old', &
|
||||||
position='append', buffered='yes')
|
position='append', buffered='yes')
|
||||||
@ -230,13 +268,13 @@ module integrals
|
|||||||
write(fname, "('errors.dat')")
|
write(fname, "('errors.dat')")
|
||||||
inquire(file=fname, exist=flag)
|
inquire(file=fname, exist=flag)
|
||||||
case default
|
case default
|
||||||
write(fname, "('errors_',i2.2,'.dat')") irun
|
write(fname, "('errors_',i2.2,'.dat')") nrun
|
||||||
flag = .false.
|
flag = .false.
|
||||||
end select
|
end select
|
||||||
|
|
||||||
! check if the file exists; if not, create a new one, otherwise move to the end
|
! check if the file exists; if not, create a new one, otherwise move to the end
|
||||||
!
|
!
|
||||||
if (flag .and. irun > 1) then
|
if (flag .and. nrun > 1) then
|
||||||
#ifdef __INTEL_COMPILER
|
#ifdef __INTEL_COMPILER
|
||||||
open(newunit=eunit, file=fname, form='formatted', status='old', &
|
open(newunit=eunit, file=fname, form='formatted', status='old', &
|
||||||
position='append', buffered='yes')
|
position='append', buffered='yes')
|
||||||
@ -290,6 +328,7 @@ module integrals
|
|||||||
subroutine finalize_integrals(status)
|
subroutine finalize_integrals(status)
|
||||||
|
|
||||||
use evolution, only : error_control
|
use evolution, only : error_control
|
||||||
|
use mpitools , only : master
|
||||||
|
|
||||||
implicit none
|
implicit none
|
||||||
|
|
||||||
@ -299,9 +338,8 @@ module integrals
|
|||||||
!
|
!
|
||||||
status = 0
|
status = 0
|
||||||
|
|
||||||
! close the integrals and statistics files
|
if (master) then
|
||||||
!
|
close(munit)
|
||||||
if (stored) then
|
|
||||||
close(funit)
|
close(funit)
|
||||||
close(sunit)
|
close(sunit)
|
||||||
if (error_control) close(eunit)
|
if (error_control) close(eunit)
|
||||||
@ -324,9 +362,12 @@ module integrals
|
|||||||
!
|
!
|
||||||
subroutine store_integrals()
|
subroutine store_integrals()
|
||||||
|
|
||||||
use blocks , only : block_meta, block_data, list_data
|
use blocks , only : block_leaf, block_data
|
||||||
|
use blocks , only : list_leaf, list_data
|
||||||
|
use blocks , only : get_mblocks, get_nleafs
|
||||||
use coordinates , only : ni => ncells, nb, ne
|
use coordinates , only : ni => ncells, nb, ne
|
||||||
use coordinates , only : advol, voli
|
use coordinates , only : advol, voli
|
||||||
|
use coordinates , only : toplev
|
||||||
use equations , only : idn, ipr, ivx, ivz, ibx, iby, ibz, ibp
|
use equations , only : idn, ipr, ivx, ivz, ibx, iby, ibz, ibp
|
||||||
use equations , only : ien, imx, imy, imz
|
use equations , only : ien, imx, imy, imz
|
||||||
use equations , only : magnetized, adiabatic_index, csnd
|
use equations , only : magnetized, adiabatic_index, csnd
|
||||||
@ -336,6 +377,8 @@ module integrals
|
|||||||
use forcing , only : einj, rinj, arms
|
use forcing , only : einj, rinj, arms
|
||||||
use helpers , only : flush_and_sync
|
use helpers , only : flush_and_sync
|
||||||
use iso_fortran_env, only : error_unit
|
use iso_fortran_env, only : error_unit
|
||||||
|
use mesh , only : changed
|
||||||
|
use mpitools , only : master, nprocs
|
||||||
#ifdef MPI
|
#ifdef MPI
|
||||||
use mpitools , only : reduce_minimum, reduce_maximum, reduce_sum
|
use mpitools , only : reduce_minimum, reduce_maximum, reduce_sum
|
||||||
#endif /* MPI */
|
#endif /* MPI */
|
||||||
@ -344,10 +387,11 @@ module integrals
|
|||||||
implicit none
|
implicit none
|
||||||
|
|
||||||
logical, save :: first = .true.
|
logical, save :: first = .true.
|
||||||
integer :: n, l, u, nk, status
|
integer :: n, l, u, nk, nl, nm, status
|
||||||
real(kind=8) :: dvol, dvolh
|
real(kind=8) :: dvol, dvolh
|
||||||
|
|
||||||
type(block_data), pointer :: pdata
|
type(block_leaf), pointer :: pleaf
|
||||||
|
type(block_data), pointer :: pdata
|
||||||
|
|
||||||
real(kind=8), dimension(:,:,:), pointer, save :: vel, mag, sqd, tmp
|
real(kind=8), dimension(:,:,:), pointer, save :: vel, mag, sqd, tmp
|
||||||
|
|
||||||
@ -358,10 +402,58 @@ module integrals
|
|||||||
real(kind=8) , parameter :: eps = epsilon(1.0d+00)
|
real(kind=8) , parameter :: eps = epsilon(1.0d+00)
|
||||||
real(kind=8) , parameter :: big = huge(1.0d+00)
|
real(kind=8) , parameter :: big = huge(1.0d+00)
|
||||||
|
|
||||||
|
integer(kind=4), dimension(toplev) :: ldist
|
||||||
|
#ifdef MPI
|
||||||
|
integer(kind=4), dimension(nprocs) :: cdist
|
||||||
|
#endif /* MPI */
|
||||||
|
|
||||||
character(len=*), parameter :: loc = 'INTEGRALS:store_integrals()'
|
character(len=*), parameter :: loc = 'INTEGRALS:store_integrals()'
|
||||||
|
|
||||||
!-------------------------------------------------------------------------------
|
!-------------------------------------------------------------------------------
|
||||||
!
|
!
|
||||||
|
! process and store the mesh statistics only on the master node
|
||||||
|
!
|
||||||
|
if (master) then
|
||||||
|
if (changed) then
|
||||||
|
|
||||||
|
! get the new numbers of meta blocks and leafs
|
||||||
|
!
|
||||||
|
nm = get_mblocks()
|
||||||
|
nl = get_nleafs()
|
||||||
|
|
||||||
|
! determine the block distribution per level and per process
|
||||||
|
!
|
||||||
|
ldist(:) = 0
|
||||||
|
#ifdef MPI
|
||||||
|
cdist(:) = 0
|
||||||
|
#endif /* MPI */
|
||||||
|
pleaf => list_leaf
|
||||||
|
do while(associated(pleaf))
|
||||||
|
l = pleaf%meta%level
|
||||||
|
ldist(l) = ldist(l) + 1
|
||||||
|
#ifdef MPI
|
||||||
|
n = pleaf%meta%process+1
|
||||||
|
cdist(n) = cdist(n) + 1
|
||||||
|
#endif /* MPI */
|
||||||
|
pleaf => pleaf%next
|
||||||
|
end do
|
||||||
|
|
||||||
|
! store the block statistics
|
||||||
|
!
|
||||||
|
write (munit,mfmt1,advance='no') step, time, fcv * nl, fef / nl, nm, nl
|
||||||
|
write (munit,mfmt2,advance='no') ldist(:)
|
||||||
|
#ifdef MPI
|
||||||
|
write (munit,mfmt3,advance='no') cdist(:)
|
||||||
|
#endif /* MPI */
|
||||||
|
write(munit,"('')")
|
||||||
|
|
||||||
|
call flush_and_sync(munit)
|
||||||
|
|
||||||
|
changed = .false.
|
||||||
|
|
||||||
|
end if ! number of blocks or leafs changed
|
||||||
|
end if ! master
|
||||||
|
|
||||||
! return if the storage interval was not reached
|
! return if the storage interval was not reached
|
||||||
!
|
!
|
||||||
if (mod(step, iintd) > 0) return
|
if (mod(step, iintd) > 0) return
|
||||||
@ -620,7 +712,7 @@ module integrals
|
|||||||
|
|
||||||
! write down the integrals and statistics to appropriate files
|
! write down the integrals and statistics to appropriate files
|
||||||
!
|
!
|
||||||
if (stored) then
|
if (master) then
|
||||||
write(funit,"(i9,13(1x,1es18.8e3))") step, time, dt, inarr(1:10), arms
|
write(funit,"(i9,13(1x,1es18.8e3))") step, time, dt, inarr(1:10), arms
|
||||||
write(sunit,"(i9,23(1x,1es18.8e3))") step, time &
|
write(sunit,"(i9,23(1x,1es18.8e3))") step, time &
|
||||||
, avarr(1), mnarr(1), mxarr(1) &
|
, avarr(1), mnarr(1), mxarr(1) &
|
||||||
|
250
sources/mesh.F90
250
sources/mesh.F90
@ -54,10 +54,6 @@ module mesh
|
|||||||
!
|
!
|
||||||
procedure(setup_problem_iface), pointer, save :: setup_problem => null()
|
procedure(setup_problem_iface), pointer, save :: setup_problem => null()
|
||||||
|
|
||||||
! the handler of the mesh statistics file
|
|
||||||
!
|
|
||||||
integer, save :: funit = 11
|
|
||||||
|
|
||||||
! the flag indicating that the block structure or distribution has changed
|
! the flag indicating that the block structure or distribution has changed
|
||||||
!
|
!
|
||||||
logical, save :: changed = .true.
|
logical, save :: changed = .true.
|
||||||
@ -70,8 +66,9 @@ module mesh
|
|||||||
|
|
||||||
public :: initialize_mesh, finalize_mesh, print_mesh
|
public :: initialize_mesh, finalize_mesh, print_mesh
|
||||||
public :: generate_mesh, update_mesh
|
public :: generate_mesh, update_mesh
|
||||||
public :: redistribute_blocks, store_mesh_stats
|
public :: redistribute_blocks
|
||||||
public :: setup_problem
|
public :: setup_problem
|
||||||
|
public :: changed
|
||||||
|
|
||||||
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
||||||
!
|
!
|
||||||
@ -92,29 +89,23 @@ module mesh
|
|||||||
!
|
!
|
||||||
! Arguments:
|
! Arguments:
|
||||||
!
|
!
|
||||||
! nrun - the restarted execution count;
|
! verbose - the verbose flag;
|
||||||
! status - the subroutine call status: 0 for success, otherwise failure;
|
! status - the subroutine call status;
|
||||||
!
|
!
|
||||||
!===============================================================================
|
!===============================================================================
|
||||||
!
|
!
|
||||||
subroutine initialize_mesh(nrun, status)
|
subroutine initialize_mesh(verbose, status)
|
||||||
|
|
||||||
use coordinates , only : nn => bcells, ni => ncells, ng => nghosts
|
use coordinates , only : ncells, nghosts
|
||||||
use coordinates , only : toplev
|
|
||||||
use equations , only : nf
|
|
||||||
use iso_fortran_env, only : error_unit
|
use iso_fortran_env, only : error_unit
|
||||||
use mpitools , only : master, nprocs
|
|
||||||
use parameters , only : get_parameter
|
|
||||||
use refinement , only : initialize_refinement
|
use refinement , only : initialize_refinement
|
||||||
|
|
||||||
implicit none
|
implicit none
|
||||||
|
|
||||||
integer, intent(in) :: nrun
|
logical, intent(in) :: verbose
|
||||||
integer, intent(out) :: status
|
integer, intent(out) :: status
|
||||||
|
|
||||||
character(len=64) :: problem = "none"
|
character(len=64) :: problem = "none"
|
||||||
character(len=64) :: fn
|
|
||||||
integer :: l, n
|
|
||||||
|
|
||||||
character(len=*), parameter :: loc = 'MESH::initialize_mesh()'
|
character(len=*), parameter :: loc = 'MESH::initialize_mesh()'
|
||||||
|
|
||||||
@ -132,53 +123,13 @@ module mesh
|
|||||||
|
|
||||||
end select
|
end select
|
||||||
|
|
||||||
! only master prepares the mesh statistics file
|
|
||||||
!
|
|
||||||
if (master) then
|
|
||||||
|
|
||||||
! generate the mesh statistics file name
|
|
||||||
!
|
|
||||||
write(fn, "('mesh_',i2.2,'.dat')") nrun
|
|
||||||
|
|
||||||
! create a new mesh statistics file
|
|
||||||
!
|
|
||||||
#ifdef __INTEL_COMPILER
|
|
||||||
open(unit = funit, file = fn, form = 'formatted', status = 'replace' &
|
|
||||||
, buffered = 'yes')
|
|
||||||
#else /* __INTEL_COMPILER */
|
|
||||||
open(unit = funit, file = fn, form = 'formatted', status = 'replace')
|
|
||||||
#endif /* __INTEL_COMPILER */
|
|
||||||
|
|
||||||
! write the mesh statistics header
|
|
||||||
!
|
|
||||||
write(funit, "('#',2(4x,a4),10x,a8,6x,a10,5x,a6,6x,a5,4x,a20)") &
|
|
||||||
'step', 'time', 'coverage', 'efficiency' &
|
|
||||||
, 'blocks', 'leafs', 'block distributions:'
|
|
||||||
|
|
||||||
! write the mesh distributions header
|
|
||||||
!
|
|
||||||
write(funit, "('#',76x,a8)", advance="no") 'level = '
|
|
||||||
do l = 1, toplev
|
|
||||||
write(funit, "(1x,i9)", advance="no") l
|
|
||||||
end do
|
|
||||||
#ifdef MPI
|
|
||||||
write(funit, "(2x,a8)", advance="no") 'nproc = '
|
|
||||||
do n = 0, nprocs - 1
|
|
||||||
write(funit, "(1x,i9)", advance="no") n
|
|
||||||
end do
|
|
||||||
#endif /* MPI */
|
|
||||||
write(funit, "('' )")
|
|
||||||
write(funit, "('#')")
|
|
||||||
|
|
||||||
end if ! master
|
|
||||||
|
|
||||||
! prepare the dimensions of the prolongation array
|
! prepare the dimensions of the prolongation array
|
||||||
!
|
!
|
||||||
pm(1:NDIMS) = 2 * (ni + ng)
|
pm(1:NDIMS) = 2 * (ncells + nghosts)
|
||||||
|
|
||||||
! initialize refinement module
|
! initialize refinement module
|
||||||
!
|
!
|
||||||
call initialize_refinement(master, status)
|
call initialize_refinement(verbose, status)
|
||||||
|
|
||||||
!-------------------------------------------------------------------------------
|
!-------------------------------------------------------------------------------
|
||||||
!
|
!
|
||||||
@ -212,8 +163,6 @@ module mesh
|
|||||||
|
|
||||||
nullify(setup_problem)
|
nullify(setup_problem)
|
||||||
|
|
||||||
if (master) close(funit)
|
|
||||||
|
|
||||||
call finalize_refinement(status)
|
call finalize_refinement(status)
|
||||||
|
|
||||||
!-------------------------------------------------------------------------------
|
!-------------------------------------------------------------------------------
|
||||||
@ -809,187 +758,6 @@ module mesh
|
|||||||
end subroutine redistribute_blocks
|
end subroutine redistribute_blocks
|
||||||
!
|
!
|
||||||
!===============================================================================
|
!===============================================================================
|
||||||
!
|
|
||||||
! subroutine STORE_MESH_STATS:
|
|
||||||
! ---------------------------
|
|
||||||
!
|
|
||||||
! Subroutine calculates and stores various mesh statistics.
|
|
||||||
!
|
|
||||||
! Arguments:
|
|
||||||
!
|
|
||||||
! step - the integration iteration step number;
|
|
||||||
! time - the physical time;
|
|
||||||
!
|
|
||||||
!===============================================================================
|
|
||||||
!
|
|
||||||
subroutine store_mesh_stats(step, time)
|
|
||||||
|
|
||||||
! import external procedures and variables
|
|
||||||
!
|
|
||||||
use blocks , only : block_meta, block_leaf
|
|
||||||
use blocks , only : list_meta, list_leaf
|
|
||||||
use blocks , only : ndims
|
|
||||||
use blocks , only : get_mblocks, get_nleafs
|
|
||||||
use coordinates , only : ddims => domain_base_dims
|
|
||||||
use coordinates , only : ni => ncells, nn => bcells
|
|
||||||
use coordinates , only : nd => nghosts_double
|
|
||||||
use coordinates , only : toplev
|
|
||||||
use helpers , only : flush_and_sync
|
|
||||||
use mpitools , only : master, nprocs
|
|
||||||
|
|
||||||
! local variables are not implicit by default
|
|
||||||
!
|
|
||||||
implicit none
|
|
||||||
|
|
||||||
! subroutine arguments
|
|
||||||
!
|
|
||||||
integer , intent(in) :: step
|
|
||||||
real(kind=8), intent(in) :: time
|
|
||||||
|
|
||||||
! local variables
|
|
||||||
!
|
|
||||||
integer(kind=4) :: l, n, nm, nl
|
|
||||||
real(kind=8) :: cv, ef, ff
|
|
||||||
|
|
||||||
! local pointers
|
|
||||||
!
|
|
||||||
type(block_meta), pointer :: pmeta
|
|
||||||
type(block_leaf), pointer :: pleaf
|
|
||||||
|
|
||||||
! local saved variables
|
|
||||||
!
|
|
||||||
logical , save :: first = .true.
|
|
||||||
real(kind=8) , save :: fcv = 1.0d+00, fef = 1.0d+00
|
|
||||||
|
|
||||||
! local arrays
|
|
||||||
!
|
|
||||||
integer(kind=4), dimension(toplev) :: ldist
|
|
||||||
#ifdef MPI
|
|
||||||
integer(kind=4), dimension(nprocs) :: cdist
|
|
||||||
#endif /* MPI */
|
|
||||||
|
|
||||||
!-------------------------------------------------------------------------------
|
|
||||||
!
|
|
||||||
! process and store the mesh statistics only on master node
|
|
||||||
!
|
|
||||||
if (master) then
|
|
||||||
|
|
||||||
! and only if the number of blocks changed
|
|
||||||
!
|
|
||||||
if (changed) then
|
|
||||||
|
|
||||||
! get the maximum number of block and level, calculated only once
|
|
||||||
!
|
|
||||||
if (first) then
|
|
||||||
|
|
||||||
! determine the number of base blocks
|
|
||||||
!
|
|
||||||
n = 0
|
|
||||||
pmeta => list_meta
|
|
||||||
|
|
||||||
do while(associated(pmeta))
|
|
||||||
|
|
||||||
if (pmeta%level == 1) n = n + 1
|
|
||||||
|
|
||||||
pmeta => pmeta%next
|
|
||||||
|
|
||||||
end do ! pmeta
|
|
||||||
|
|
||||||
! calculate the maximum number of blocks
|
|
||||||
!
|
|
||||||
n = n * 2**(ndims*(toplev - 1))
|
|
||||||
|
|
||||||
! calculate the coverage and efficiency factors
|
|
||||||
!
|
|
||||||
ff = 2**(toplev - 1)
|
|
||||||
fcv = 1.0d+00 / n
|
|
||||||
fef = 1.0d+00 * (ddims(1) * ni * ff + nd) / nn
|
|
||||||
fef = fef * (ddims(2) * ni * ff + nd) / nn
|
|
||||||
#if NDIMS == 3
|
|
||||||
fef = fef * (ddims(3) * ni * ff + nd) / nn
|
|
||||||
#endif /* NDIMS == 3 */
|
|
||||||
|
|
||||||
! reset the first execution flag
|
|
||||||
!
|
|
||||||
first = .false.
|
|
||||||
|
|
||||||
end if ! first
|
|
||||||
|
|
||||||
! get the new numbers of meta blocks and leafs
|
|
||||||
!
|
|
||||||
nm = get_mblocks()
|
|
||||||
nl = get_nleafs()
|
|
||||||
|
|
||||||
! calculate the coverage (the number of leafs divided by the maximum number
|
|
||||||
! of blocks) and the efficiency (the total number of cells in a corresponding
|
|
||||||
! uniform mesh divided by the number of cells for adaptive mesh case)
|
|
||||||
!
|
|
||||||
cv = fcv * nl
|
|
||||||
ef = fef / nl
|
|
||||||
|
|
||||||
! get the block level and block process distributions
|
|
||||||
!
|
|
||||||
ldist(:) = 0
|
|
||||||
#ifdef MPI
|
|
||||||
cdist(:) = 0
|
|
||||||
#endif /* MPI */
|
|
||||||
|
|
||||||
pleaf => list_leaf
|
|
||||||
|
|
||||||
do while(associated(pleaf))
|
|
||||||
|
|
||||||
pmeta => pleaf%meta
|
|
||||||
|
|
||||||
ldist(pmeta%level) = ldist(pmeta%level) + 1
|
|
||||||
#ifdef MPI
|
|
||||||
cdist(pmeta%process+1) = cdist(pmeta%process+1) + 1
|
|
||||||
#endif /* MPI */
|
|
||||||
|
|
||||||
pleaf => pleaf%next
|
|
||||||
|
|
||||||
end do ! pleaf
|
|
||||||
|
|
||||||
! store the block statistics
|
|
||||||
!
|
|
||||||
write(funit, "(i9,3es14.6e3,2(2x,i9))", advance="no") &
|
|
||||||
step, time, cv, ef, nm, nl
|
|
||||||
|
|
||||||
! store the block level distribution
|
|
||||||
!
|
|
||||||
write(funit,"(12x)", advance="no")
|
|
||||||
do l = 1, toplev
|
|
||||||
write(funit,"(1x,i9)", advance="no") ldist(l)
|
|
||||||
end do ! l = 1, toplev
|
|
||||||
|
|
||||||
#ifdef MPI
|
|
||||||
! store the process level distribution
|
|
||||||
!
|
|
||||||
write(funit,"(10x)", advance="no")
|
|
||||||
do l = 1, nprocs
|
|
||||||
write(funit,"(1x,i9)", advance="no") cdist(l)
|
|
||||||
end do ! l = 1, nprocs
|
|
||||||
#endif /* MPI */
|
|
||||||
|
|
||||||
! append the new line character
|
|
||||||
!
|
|
||||||
write(funit,"('')")
|
|
||||||
|
|
||||||
call flush_and_sync(funit)
|
|
||||||
|
|
||||||
end if ! number of blocks or leafs changed
|
|
||||||
|
|
||||||
end if ! master
|
|
||||||
|
|
||||||
! reset the flag indicating the block structure has changed
|
|
||||||
!
|
|
||||||
changed = .false.
|
|
||||||
|
|
||||||
|
|
||||||
!-------------------------------------------------------------------------------
|
|
||||||
!
|
|
||||||
end subroutine store_mesh_stats
|
|
||||||
!
|
|
||||||
!===============================================================================
|
|
||||||
!!
|
!!
|
||||||
!!*** PRIVATE SUBROUTINES ****************************************************
|
!!*** PRIVATE SUBROUTINES ****************************************************
|
||||||
!!
|
!!
|
||||||
|
@ -178,7 +178,7 @@ module problem
|
|||||||
return
|
return
|
||||||
end if
|
end if
|
||||||
|
|
||||||
call initialize_mesh(nrun, status)
|
call initialize_mesh(verbose, status)
|
||||||
if (status /= 0) then
|
if (status /= 0) then
|
||||||
write(error_unit,"('[',a,']: ',a)") trim(loc), &
|
write(error_unit,"('[',a,']: ',a)") trim(loc), &
|
||||||
"Could not initialize mesh!"
|
"Could not initialize mesh!"
|
||||||
|
Loading…
x
Reference in New Issue
Block a user