Particle release scheduling.
◆ advance()
subroutine particlereleaseschedulemodule::advance |
( |
class(particlereleasescheduletype), intent(inout) |
this, |
|
|
character(len=linelength), dimension(:), intent(in), optional |
lines |
|
) |
| |
|
private |
This involves several tasks: first, advance the time selection. Then, if period-block release setting lines are provided, reinitialize the time step selection for the given period. Finally, refresh the schedule array, deduplicating any times closer than the set tolerance.
This routine is idempotent.
Definition at line 108 of file ParticleReleaseSchedule.f90.
110 class(ParticleReleaseScheduleType),
intent(inout) :: this
111 character(len=LINELENGTH),
intent(in),
optional :: lines(:)
112 integer(I4B) :: it, i
113 real(DP) :: tprevious
117 call this%time_select%advance()
121 if (
present(lines))
then
122 call this%step_select%init()
123 do i = 1,
size(lines)
124 call this%step_select%read(lines(i))
129 if (
allocated(this%times))
deallocate (this%times)
130 allocate (this%times(0))
139 call this%schedule(trelease)
145 if (this%time_select%any())
then
146 do it = this%time_select%selection(1), this%time_select%selection(2)
147 trelease = this%time_select%times(it)
148 if (tprevious >= dzero .and. is_close( &
151 atol=this%tolerance)) cycle
153 call this%schedule(trelease)
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
real(dp), pointer, public totimc
simulation time at start of time step
integer(i4b), pointer, public kstp
current time step number
◆ any()
Note: be sure to call advance() before calling this function, or the result may still be associated with a prior time step.
Definition at line 164 of file ParticleReleaseSchedule.f90.
165 class(ParticleReleaseScheduleType) :: this
◆ count()
Note: be sure to call advance() before calling this function, or the result may still be associated with a prior time step.
Definition at line 174 of file ParticleReleaseSchedule.f90.
175 class(ParticleReleaseScheduleType) :: this
◆ create_release_schedule()
type(particlereleasescheduletype) function, pointer, public particlereleaseschedulemodule::create_release_schedule |
( |
real(dp), intent(in) |
tolerance | ) |
|
- Parameters
-
[in] | tolerance | coincident release time tolerance |
- Returns
- schedule pointer
Definition at line 46 of file ParticleReleaseSchedule.f90.
47 real(DP),
intent(in) :: tolerance
48 type(ParticleReleaseScheduleType),
pointer :: schedule_
51 allocate (schedule_%times(0))
52 allocate (schedule_%time_select)
53 allocate (schedule_%step_select)
54 call schedule_%time_select%init()
55 call schedule_%step_select%init()
56 schedule_%tolerance = tolerance
◆ destroy()
- Parameters
-
[in,out] | this | this instance |
Definition at line 60 of file ParticleReleaseSchedule.f90.
61 class(ParticleReleaseScheduleType),
intent(inout) :: this
63 deallocate (this%times)
64 call this%time_select%deallocate()
65 call this%step_select%deallocate()
66 deallocate (this%time_select)
67 deallocate (this%step_select)
◆ log()
subroutine particlereleaseschedulemodule::log |
( |
class(particlereleasescheduletype), intent(inout) |
this, |
|
|
integer(i4b), intent(in) |
iout |
|
) |
| |
|
private |
- Parameters
-
[in,out] | this | this instance |
[in] | iout | output unit |
Definition at line 71 of file ParticleReleaseSchedule.f90.
72 class(ParticleReleaseScheduleType),
intent(inout) :: this
73 integer(I4B),
intent(in) :: iout
74 character(len=*),
parameter :: fmt = &
75 &
"(6x,A,': ',50(G0,' '))"
78 write (iout, fmt)
'RELEASE SCHEDULE', this%times
80 write (iout,
"(1x,a,1x,a)")
'NO RELEASES SCHEDULED'
◆ schedule()
subroutine particlereleaseschedulemodule::schedule |
( |
class(particlereleasescheduletype), intent(inout) |
this, |
|
|
real(dp), intent(in) |
trelease |
|
) |
| |
|
private |
To schedule multiple release times at once, expand and populate the time selection object by hand. DO NOT attempt to manipulate the times array; this is a read-only property which the schedule maintains.
Definition at line 91 of file ParticleReleaseSchedule.f90.
92 class(ParticleReleaseScheduleType),
intent(inout) :: this
93 real(DP),
intent(in) :: trelease
94 call expandarray(this%times)
95 this%times(
size(this%times)) = trelease