2010-10-06 23:03:47 -03:00
|
|
|
!!******************************************************************************
|
2008-11-04 13:08:01 -06:00
|
|
|
!!
|
|
|
|
!! program: Godunov-AMR
|
|
|
|
!!
|
2011-02-27 22:45:54 -03:00
|
|
|
!! Copyright (C) 2008-2011 Grzegorz Kowal <grzegorz@gkowal.info>
|
2008-11-04 13:08:01 -06:00
|
|
|
!!
|
2010-10-06 23:03:47 -03:00
|
|
|
!!******************************************************************************
|
2008-11-04 13:08:01 -06:00
|
|
|
!!
|
|
|
|
!! This file is part of Godunov-AMR.
|
|
|
|
!!
|
|
|
|
!! Godunov-AMR is free software; you can redistribute it and/or modify
|
|
|
|
!! it under the terms of the GNU General Public License as published by
|
|
|
|
!! the Free Software Foundation; either version 3 of the License, or
|
|
|
|
!! (at your option) any later version.
|
|
|
|
!!
|
|
|
|
!! Godunov-AMR is distributed in the hope that it will be useful,
|
|
|
|
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
!! GNU General Public License for more details.
|
|
|
|
!!
|
|
|
|
!! You should have received a copy of the GNU General Public License
|
|
|
|
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
!!
|
2010-10-06 23:03:47 -03:00
|
|
|
!!******************************************************************************
|
2008-11-04 13:08:01 -06:00
|
|
|
!!
|
|
|
|
!
|
|
|
|
program godunov
|
|
|
|
|
2008-11-04 21:00:50 -06:00
|
|
|
! modules
|
|
|
|
!
|
2011-02-28 16:04:38 -03:00
|
|
|
use blocks , only : nleafs
|
2011-04-14 00:19:45 +02:00
|
|
|
use blocks , only : init_blocks
|
2011-04-11 16:27:08 +02:00
|
|
|
use config , only : read_config, nres, nmax, tmax, dtini, dtout, ftype, cfl
|
2010-12-08 18:11:15 -02:00
|
|
|
#ifdef FORCE
|
|
|
|
use config , only : fdt
|
|
|
|
#endif /* FORCE */
|
2010-12-01 10:43:01 -02:00
|
|
|
use evolution, only : evolve, update_maximum_speed, n, t, dt, dtn
|
2010-12-08 12:21:41 -02:00
|
|
|
#ifdef FORCE
|
|
|
|
use forcing , only : init_forcing, clear_forcing
|
|
|
|
#endif /* FORCE */
|
2011-03-02 15:52:33 -03:00
|
|
|
use integrals, only : init_integrals, clear_integrals, store_integrals
|
2011-04-11 16:27:08 +02:00
|
|
|
use io , only : write_data, restart_job
|
2008-12-07 18:57:08 -06:00
|
|
|
use mesh , only : init_mesh, clear_mesh
|
2008-12-22 15:09:05 -06:00
|
|
|
use mpitools , only : ncpu, ncpus, init_mpi, clear_mpi, is_master
|
2010-12-08 11:37:56 -02:00
|
|
|
use random , only : init_generator
|
2010-10-06 23:11:22 -03:00
|
|
|
use timer , only : init_timers, start_timer, stop_timer, get_timer &
|
|
|
|
, get_timer_total
|
2008-11-04 13:08:01 -06:00
|
|
|
!
|
2008-12-08 21:11:17 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-04 13:08:01 -06:00
|
|
|
!
|
|
|
|
! local variables
|
|
|
|
!
|
2008-12-07 14:06:04 -06:00
|
|
|
character(len=60) :: fmt
|
2008-12-16 22:34:54 -06:00
|
|
|
integer :: no, ed, eh, em, es, ec
|
2010-10-06 23:11:22 -03:00
|
|
|
real :: tall, tbeg, tcur, per
|
2008-11-04 13:08:01 -06:00
|
|
|
!
|
2008-12-08 21:11:17 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-04 21:00:50 -06:00
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
! initialize MPI
|
|
|
|
!
|
2010-12-08 11:46:09 -02:00
|
|
|
call init_mpi()
|
2008-12-22 14:57:31 -06:00
|
|
|
|
2008-11-04 21:00:50 -06:00
|
|
|
! print info message
|
2008-11-04 13:08:01 -06:00
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
if (is_master()) then
|
|
|
|
write (*,"(1x,78('-'))")
|
|
|
|
write (*,"(1x,18('='),4x,a,4x,19('='))") ' Godunov-AMR algorithm '
|
2011-03-21 18:24:55 -03:00
|
|
|
write (*,"(1x,16('='),4x,a,4x,16('='))") &
|
|
|
|
'Copyright (C) 2008-2011 Grzegorz Kowal'
|
2008-12-22 15:09:05 -06:00
|
|
|
#ifdef MPI
|
2010-10-06 23:11:22 -03:00
|
|
|
write (*,"(1x,18('='),4x,a,i5,a,4x,19('='))") 'MPI enabled with ', ncpus &
|
|
|
|
, ' processors'
|
2008-12-22 15:09:05 -06:00
|
|
|
#endif /* MPI */
|
2008-12-22 14:57:31 -06:00
|
|
|
write (*,"(1x,78('-'))")
|
|
|
|
write (*,*)
|
2011-03-21 18:24:55 -03:00
|
|
|
write (*,"(1x,a)" ) "Physics:"
|
|
|
|
write (*,"(4x,a,1x,a)" ) "equations =", &
|
|
|
|
#ifdef HYDRO
|
|
|
|
"HD"
|
|
|
|
#endif /* HYDRO */
|
|
|
|
#ifdef MHD
|
|
|
|
"MHD"
|
|
|
|
#endif /* MHD */
|
|
|
|
write (*,"(4x,a,1x,a)" ) "equation of state =", &
|
|
|
|
#ifdef ADI
|
2011-03-22 15:36:59 -03:00
|
|
|
"adiabatic"
|
2011-03-21 18:24:55 -03:00
|
|
|
#endif /* ADI */
|
|
|
|
#ifdef ISO
|
|
|
|
"isothermal"
|
|
|
|
#endif /* ISO */
|
|
|
|
|
|
|
|
write (*,*)
|
2010-10-06 23:11:22 -03:00
|
|
|
end if
|
2008-11-04 13:08:01 -06:00
|
|
|
|
2008-11-05 22:16:24 -06:00
|
|
|
! read configuration file
|
|
|
|
!
|
2010-12-08 11:46:09 -02:00
|
|
|
call read_config()
|
2008-11-05 22:16:24 -06:00
|
|
|
|
2011-04-14 00:19:45 +02:00
|
|
|
! reset number of iterations and time, etc.
|
|
|
|
!
|
|
|
|
n = 0
|
|
|
|
t = 0.0
|
|
|
|
dt = cfl * dtini
|
|
|
|
dtn = dtini
|
|
|
|
no = 0
|
|
|
|
tbeg = 0.0
|
|
|
|
|
2008-12-07 14:06:04 -06:00
|
|
|
! initialize timers
|
|
|
|
!
|
2010-12-08 11:46:09 -02:00
|
|
|
call init_timers()
|
2008-12-07 14:06:04 -06:00
|
|
|
|
2010-12-08 12:21:41 -02:00
|
|
|
! initialize random number generator
|
2010-12-08 11:37:56 -02:00
|
|
|
!
|
|
|
|
call init_generator()
|
|
|
|
|
2011-04-14 00:19:45 +02:00
|
|
|
! initialize block module
|
2008-12-07 15:14:18 -06:00
|
|
|
!
|
2011-04-14 00:19:45 +02:00
|
|
|
call start_timer(1)
|
|
|
|
call init_blocks()
|
|
|
|
call stop_timer(1)
|
2011-03-02 15:04:16 -03:00
|
|
|
|
2011-04-11 16:27:08 +02:00
|
|
|
! check if we initiate new problem or restart previous job
|
|
|
|
!
|
|
|
|
if (nres .lt. 0) then
|
|
|
|
|
|
|
|
! initialize our adaptive mesh, refine that mesh to the desired level
|
|
|
|
! according to the initialized problem
|
|
|
|
!
|
|
|
|
call start_timer(1)
|
|
|
|
call init_mesh()
|
|
|
|
call stop_timer(1)
|
|
|
|
|
|
|
|
else
|
|
|
|
|
|
|
|
! reconstruct the meta and data block structures from a given restart file
|
|
|
|
!
|
|
|
|
call start_timer(1)
|
|
|
|
call restart_job(nres, ncpu)
|
|
|
|
call stop_timer(1)
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
2010-12-08 18:11:15 -02:00
|
|
|
#ifdef FORCE
|
|
|
|
! if the forcing time step is larger than the initial time step decrease it
|
|
|
|
!
|
|
|
|
fdt = min(fdt, 0.1d0 * dtini)
|
|
|
|
|
|
|
|
! round the initial time step to the integer number of forcing time steps
|
|
|
|
!
|
|
|
|
dt = fdt * floor(dt / fdt)
|
2008-12-07 15:14:18 -06:00
|
|
|
|
2011-03-02 15:19:01 -03:00
|
|
|
! initialize forcing module
|
|
|
|
!
|
|
|
|
call init_forcing()
|
|
|
|
|
2010-12-08 12:21:41 -02:00
|
|
|
#endif /* FORCE */
|
2011-03-02 15:52:33 -03:00
|
|
|
! initialize the integrals module
|
|
|
|
!
|
|
|
|
call init_integrals()
|
|
|
|
|
2010-12-01 10:43:01 -02:00
|
|
|
! update the maximum speed in the system
|
|
|
|
!
|
|
|
|
call start_timer(2)
|
|
|
|
call update_maximum_speed()
|
|
|
|
call stop_timer(2)
|
|
|
|
|
2008-11-29 22:22:19 -06:00
|
|
|
! write down the initial state
|
|
|
|
!
|
2008-12-07 18:57:08 -06:00
|
|
|
call start_timer(3)
|
2008-12-22 14:57:31 -06:00
|
|
|
call write_data(ftype, no, ncpu)
|
2008-12-07 18:57:08 -06:00
|
|
|
call stop_timer(3)
|
2008-11-29 22:22:19 -06:00
|
|
|
|
2008-12-07 15:14:18 -06:00
|
|
|
! print information
|
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
if (is_master()) then
|
2011-03-10 15:22:10 -03:00
|
|
|
write(*,"(1x,a)" ) "Evolving the system:"
|
2011-02-28 16:04:38 -03:00
|
|
|
write(*,'(4x,a4,5x,a4,11x,a2,12x,a6,7x,a3)') 'step', 'time', 'dt' &
|
|
|
|
, 'blocks', 'ETA'
|
|
|
|
write(*,'(i8,2(1x,1pe14.6),2x,i8,2x,1i4.1,"d",1i2.2,"h",1i2.2,"m"' // &
|
|
|
|
',1i2.2,"s",a1,$)') n, t, dt, nleafs, ed, eh, em, es, char(13)
|
2010-10-06 23:11:22 -03:00
|
|
|
end if
|
2008-12-07 15:14:18 -06:00
|
|
|
|
|
|
|
! main loop
|
|
|
|
!
|
|
|
|
do while((n .lt. nmax) .and. (t .le. tmax))
|
|
|
|
|
2008-12-09 14:51:33 -06:00
|
|
|
! compute new time step
|
|
|
|
!
|
2010-12-08 18:11:15 -02:00
|
|
|
dt = min(cfl * dtn, 2.0d0 * dt)
|
|
|
|
|
|
|
|
#ifdef FORCE
|
|
|
|
! limit the time step to the integer number of forcing time steps
|
|
|
|
!
|
|
|
|
dt = fdt * floor(dt / fdt)
|
|
|
|
#endif /* FORCE */
|
2008-12-09 14:51:33 -06:00
|
|
|
|
2008-12-07 15:14:18 -06:00
|
|
|
! advance the iteration number and time
|
|
|
|
!
|
|
|
|
t = t + dt
|
|
|
|
n = n + 1
|
|
|
|
|
2008-12-07 18:57:08 -06:00
|
|
|
! performe one step evolution
|
|
|
|
!
|
|
|
|
call start_timer(2)
|
2010-12-08 11:46:09 -02:00
|
|
|
call evolve()
|
2008-12-07 18:57:08 -06:00
|
|
|
call stop_timer(2)
|
|
|
|
|
2011-03-02 15:52:33 -03:00
|
|
|
! store integrals
|
|
|
|
!
|
|
|
|
call store_integrals()
|
|
|
|
|
2008-12-09 12:42:10 -06:00
|
|
|
! store data
|
|
|
|
!
|
|
|
|
call start_timer(3)
|
|
|
|
if (dtout .gt. 0.0 .and. no .lt. (int(t/dtout))) then
|
|
|
|
no = no + 1
|
2008-12-22 14:57:31 -06:00
|
|
|
call write_data(ftype, no, ncpu)
|
2010-10-06 23:11:22 -03:00
|
|
|
end if
|
2008-12-09 12:42:10 -06:00
|
|
|
call stop_timer(3)
|
|
|
|
|
2008-12-16 22:34:54 -06:00
|
|
|
! get current time in seconds
|
|
|
|
!
|
|
|
|
tcur = get_timer_total()
|
|
|
|
ec = int((tmax - t)/(t - tbeg)*tcur, kind=4)
|
|
|
|
es = max(0, int(mod(ec,60)))
|
|
|
|
em = int(mod(ec/60,60))
|
|
|
|
eh = int(ec/3600)
|
|
|
|
ed = int(eh/24)
|
|
|
|
eh = int(mod(eh, 24))
|
|
|
|
ed = min(9999,ed)
|
|
|
|
|
2008-12-07 15:14:18 -06:00
|
|
|
! print progress information
|
|
|
|
!
|
2011-02-27 22:39:56 -03:00
|
|
|
if (is_master()) &
|
2011-02-28 16:04:38 -03:00
|
|
|
write(*,'(i8,2(1x,1pe14.6),2x,i8,2x,1i4.1,"d",1i2.2,"h",1i2.2,"m"' // &
|
|
|
|
',1i2.2,"s",a1,$)') n, t, dt, nleafs, ed, eh, em, es, char(13)
|
2008-12-07 15:14:18 -06:00
|
|
|
end do
|
|
|
|
|
2011-02-27 22:39:56 -03:00
|
|
|
! add one empty line
|
|
|
|
!
|
|
|
|
if (is_master()) write(*,*)
|
|
|
|
|
2008-12-08 21:14:12 -06:00
|
|
|
! write down the final state
|
|
|
|
!
|
|
|
|
call start_timer(3)
|
2008-12-09 12:42:10 -06:00
|
|
|
no = no + 1
|
2008-12-22 14:57:31 -06:00
|
|
|
call write_data(ftype, no, ncpu)
|
2008-12-08 21:14:12 -06:00
|
|
|
call stop_timer(3)
|
|
|
|
|
2010-12-08 12:21:41 -02:00
|
|
|
#ifdef FORCE
|
|
|
|
! finalize forcing module
|
|
|
|
!
|
|
|
|
call clear_forcing()
|
2011-03-02 15:04:16 -03:00
|
|
|
|
2010-12-08 12:21:41 -02:00
|
|
|
#endif /* FORCE */
|
2011-03-02 15:52:33 -03:00
|
|
|
! clear up the integrals module
|
|
|
|
!
|
|
|
|
call clear_integrals()
|
|
|
|
|
2011-03-02 15:04:16 -03:00
|
|
|
! deallocate and reset mesh
|
|
|
|
!
|
|
|
|
call start_timer(1)
|
|
|
|
call clear_mesh()
|
|
|
|
call stop_timer(1)
|
2010-12-08 12:21:41 -02:00
|
|
|
|
2008-12-07 14:06:04 -06:00
|
|
|
! get total time
|
|
|
|
!
|
|
|
|
tall = get_timer_total()
|
2010-10-06 23:11:22 -03:00
|
|
|
per = 100.0 / tall
|
2008-12-07 14:06:04 -06:00
|
|
|
|
|
|
|
! print info about execution times
|
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
if (is_master()) then
|
2011-03-10 15:22:10 -03:00
|
|
|
write(fmt,"(a,i2,a)") "(a27,1f", max(1, nint(alog10(tall))) + 5 &
|
|
|
|
, ".3,' secs = ',f6.2,' %')"
|
2008-12-07 14:06:04 -06:00
|
|
|
|
2008-12-22 14:57:31 -06:00
|
|
|
write (*,*)
|
2011-03-10 15:22:10 -03:00
|
|
|
write(*,"(1x,a)") "Job timings:"
|
|
|
|
write (*,fmt) "Initialization : ", get_timer(1), per * get_timer(1)
|
|
|
|
write (*,fmt) "Evolution : ", get_timer(2), per * get_timer(2)
|
|
|
|
write (*,fmt) "Data output : ", get_timer(3), per * get_timer(3)
|
|
|
|
write (*,fmt) "Boundary update : ", get_timer(4), per * get_timer(4)
|
|
|
|
write (*,fmt) "Mesh update : ", get_timer(5), per * get_timer(5)
|
2011-03-01 18:58:18 -03:00
|
|
|
#ifdef FORCE
|
2011-03-10 15:22:10 -03:00
|
|
|
write (*,fmt) "External forcing : ", get_timer(10), per * get_timer(10)
|
|
|
|
write (*,fmt) " - initialization : ", get_timer(11), per * get_timer(11)
|
|
|
|
write (*,fmt) " - evolution : ", get_timer(12), per * get_timer(12)
|
|
|
|
write (*,fmt) " - real to fourier : ", get_timer(13), per * get_timer(13)
|
|
|
|
write (*,fmt) " - fourier to real : ", get_timer(14), per * get_timer(14)
|
2011-03-01 18:58:18 -03:00
|
|
|
#endif /* FORCE */
|
2011-03-10 15:22:10 -03:00
|
|
|
write (*,fmt) "EXECUTION TIME : ", tall , 100.0
|
|
|
|
write (*,*)
|
2010-10-06 23:11:22 -03:00
|
|
|
end if
|
2008-12-22 14:57:31 -06:00
|
|
|
|
|
|
|
! close access to the MPI
|
|
|
|
!
|
2010-12-08 11:46:09 -02:00
|
|
|
call clear_mpi()
|
2008-12-07 14:06:04 -06:00
|
|
|
|
2008-12-08 21:11:17 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-11-04 13:08:01 -06:00
|
|
|
!
|
|
|
|
end program
|