-
Notifications
You must be signed in to change notification settings - Fork 2
/
timestepmod.f90
127 lines (112 loc) · 4.38 KB
/
timestepmod.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
! © 2023. Triad National Security, LLC. All rights reserved.
! This program was produced under U.S. Government contract 89233218CNA000001 for Los Alamos National
! Laboratory (LANL), which is operated by Triad National Security, LLC for the U.S. Department of
! Energy/National Nuclear Security Administration. All rights in the program are reserved by Triad
! National Security, LLC, and the U.S. Department of Energy/National Nuclear Security Administration.
! The Government is granted for itself and others acting on its behalf a nonexclusive, paid-up,
! irrevocable worldwide license in this material to reproduce, prepare. derivative works, distribute
! copies to the public, perform publicly and display publicly, and to permit others to do so.
!This file is part of SuperNu. SuperNu is released under the terms of the GNU GPLv3, see COPYING.
!Copyright (c) 2013-2022 Ryan T. Wollaeger and Daniel R. van Rossum. All rights reserved.
module timestepmod
implicit none
character(4) :: tsp_gridtype = 'none' !inputpar
integer :: tsp_nt = 0 !total # of time steps
integer :: tsp_itrestart = 0 !restart time step # (at beginning of time step)
integer :: tsp_it !current time step number
real*8 :: tsp_t,tsp_t1
real*8,allocatable :: tsp_tarr(:) !time step array
real*8,allocatable :: tsp_tpreset(:) !store preset time steps from input.tsp_time
real*8 :: tsp_tcenter,tsp_tfirst,tsp_tlast
real*8 :: tsp_dt,tsp_dtinv
private read_timestep_preset
save
contains
subroutine timestepmod_init
!----------------------------!{{{
use physconstmod
!***********************************************************************
! set the timestep constants
!***********************************************************************
!-- read timestep configuration from file
if(tsp_gridtype=='read') then
call read_timestep_preset(tsp_nt)
tsp_t = tsp_tpreset(1)
else
!-- configured by input parameters
tsp_t = tsp_tfirst
endif
!-- alloc
allocate(tsp_tarr(tsp_nt+1))
tsp_tarr(1) = tsp_t
!!}}}
end subroutine timestepmod_init
subroutine timestep_update
implicit none!{{{
!***********************************************************************
! update the timestep variables
!***********************************************************************
real*8 :: help
!
!-- preset time step sizes
select case(tsp_gridtype)
case('read')
if(.not.allocated(tsp_tpreset)) stop 'timestep_update: tpreset not allocated'
tsp_t = tsp_tpreset(tsp_it)
tsp_t1 = tsp_tpreset(tsp_it+1)
tsp_dt = tsp_t1 - tsp_t
case('lin ')
!-- linear time grid
tsp_dt = (tsp_tlast - tsp_tfirst)/tsp_nt
tsp_t1 = tsp_tfirst + tsp_it*tsp_dt !beginning of the time step
tsp_t = tsp_tfirst + (tsp_it-1)*tsp_dt !beginning of the time step
case('expo')
!-- exponential time grid
help = log(tsp_tlast/tsp_tfirst)/tsp_nt
tsp_t1 = tsp_tfirst*exp(tsp_it*help) !beginning of the time step
tsp_t = tsp_tfirst*exp((tsp_it-1)*help) !beginning of the time step
tsp_dt = tsp_t1 - tsp_t
case default
stop 'timestep_update: invalid tsp_gridtype'
end select
!-- append in time array
tsp_tarr(tsp_it+1) = tsp_t1
tsp_tcenter = tsp_t + .5*tsp_dt
tsp_dtinv = 1d0/tsp_dt
!}}}
end subroutine timestep_update
subroutine read_timestep_preset(nt)
implicit none!{{{
integer,intent(out) :: nt
!***********************************************************************
! read preset timestep values
!***********************************************************************
integer :: istat
character :: dmy
!
open(4,file='input.tsp_time',status='old',iostat=istat)
if(istat/=0) stop 'rd_tsp_preset: file missing: input.tsp_time'
!
!-- read header
read(4,*,iostat=istat) dmy, nt
if(istat/=0) stop 'rd_tsp_preset: file header error: input.tsp_time'
!
!-- allocate
if(nt<1) stop 'rd_tsp_preset: header nt < 1'
allocate(tsp_tpreset(nt+1))
!
!-- read lines
read(4,*,iostat=istat) tsp_tpreset
if(istat/=0) stop 'rd_tsp_preset: file error: input.tsp_time'
!
!-- make sure no remaining data
read(4,*,iostat=istat) !-- last time step value is only for tsp_dt purpose
if(istat==0) stop 'rd_tsp_preset: file body too long'
!
close(4)
!
write(6,*) 'rd_tsp_preset: custom timesteps read successfully',nt
!!}}}
end subroutine read_timestep_preset
end module timestepmod
! vim: fdm=marker