40 double precision :: sts_dtpar=0.5d0
43 double precision,
parameter :: nu_sts = 0.5d0
45 integer :: sts_ncycles=1000
46 integer :: sts_method = 1
52 logical :: fix_conserve_at_step = .true.
53 logical :: sts_initialized = .false.
58 subroutine subr1(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
60 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
61 double precision,
intent(in) :: x(ixi^s,1:
ndim)
62 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
63 double precision,
intent(in) :: my_dt
64 logical,
intent(in) :: fix_conserve_at_step
68 function subr2(w,ixG^L,ix^L,dx^D,x)
result(dtnew)
70 integer,
intent(in) :: ixg^
l, ix^
l
71 double precision,
intent(in) ::
dx^
d, x(ixg^s,1:
ndim)
72 double precision,
intent(in) :: w(ixg^s,1:nw)
73 double precision :: dtnew
77 subroutine subr_e(w, x, ixI^L, ixO^L, step)
80 integer,
intent(in) :: ixi^
l,ixo^
l
81 double precision,
intent(inout) :: w(ixi^s,1:nw)
82 double precision,
intent(in) :: x(ixi^s,1:
ndim)
83 integer,
intent(in) :: step
87 subroutine subr5(ixI^L, ixO^L, w, x)
89 integer,
intent(in) :: ixi^
l, ixo^
l
90 double precision,
intent(in) :: x(ixi^s,1:
ndim)
91 double precision,
intent(inout) :: w(ixi^s,1:nw)
98 double precision,
intent(in) :: dt
102 function subr4(dt,dtnew,dt_modified)
result(s)
103 double precision,
intent(in) :: dtnew
104 double precision,
intent(inout) :: dt
105 logical,
intent(inout) :: dt_modified
113 double precision :: dt_expl
117 integer,
dimension(-1:1^D&) :: type_send_srl_sts_1, type_recv_srl_sts_1
118 integer,
dimension(-1:1^D&) :: type_send_r_sts_1
119 integer,
dimension( 0:3^D&) :: type_recv_r_sts_1
120 integer,
dimension( 0:3^D&) :: type_recv_p_sts_1, type_send_p_sts_1
122 integer,
dimension(-1:1^D&) :: type_send_srl_sts_2, type_recv_srl_sts_2
123 integer,
dimension(-1:1^D&) :: type_send_r_sts_2
124 integer,
dimension( 0:3^D&) :: type_recv_r_sts_2
125 integer,
dimension( 0:3^D&) :: type_recv_p_sts_2, type_send_p_sts_2
133 logical :: types_initialized
134 logical :: evolve_magnetic_field
135 procedure(subr1),
pointer,
nopass :: sts_set_sources
136 procedure(subr2),
pointer,
nopass :: sts_getdt
137 procedure(subr5),
pointer,
nopass :: sts_before_first_cycle, sts_after_last_cycle
138 procedure(subr_e),
pointer,
nopass :: sts_handle_errors
139 type(sts_term),
pointer :: next
147 procedure(subr4),
pointer :: sts_get_ncycles
155 if(.not. sts_initialized)
then
158 sts_dtpar=sts_dtpar/dble(
ndim)
159 sts_initialized = .true.
160 if(sts_method .eq. 1)
then
162 sts_get_ncycles => sts_get_ncycles1
163 else if(sts_method .eq. 2)
then
165 sts_get_ncycles => sts_get_ncycles2
167 call mpistop(
"Unknown sts method")
175 if (sts_initialized)
then
183 subroutine sts_params_read(files)
185 character(len=*),
intent(in) :: files(:)
188 namelist /sts_list/ sts_dtpar,sts_ncycles,sts_method,
sourcetype_sts
190 do n = 1,
size(files)
191 open(
unitpar, file=trim(files(n)), status=
"old")
192 read(
unitpar, sts_list,
end=111)
196 end subroutine sts_params_read
205 subroutine add_sts_method(sts_getdt, sts_set_sources, startVar, nflux, startwbc, nwbc, evolve_B)
209 integer,
intent(in) :: startvar, nflux, startwbc, nwbc
210 logical,
intent(in) :: evolve_b
214 subroutine sts_set_sources(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
217 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
218 double precision,
intent(in) :: x(ixi^s,1:
ndim)
219 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
220 double precision,
intent(in) :: my_dt
221 logical,
intent(in) :: fix_conserve_at_step
222 end subroutine sts_set_sources
224 function sts_getdt(w,ixG^L,ix^L,dx^D,x)
result(dtnew)
226 integer,
intent(in) :: ixg^
l, ix^
l
227 double precision,
intent(in) ::
dx^
d, x(ixg^s,1:
ndim)
228 double precision,
intent(in) :: w(ixg^s,1:nw)
229 double precision :: dtnew
230 end function sts_getdt
234 type(sts_term),
pointer :: temp
237 temp%sts_getdt => sts_getdt
238 temp%sts_set_sources => sts_set_sources
239 temp%sts_before_first_cycle => null()
240 temp%sts_after_last_cycle => null()
241 temp%sts_handle_errors => null()
242 temp%startVar = startvar
243 temp%endVar= startvar+nflux-1
245 temp%startwbc = startwbc
247 temp%types_initialized = .false.
248 temp%evolve_magnetic_field=evolve_b
260 subroutine sts_before_first_cycle(ixI^L, ixO^L, w, x)
262 integer,
intent(in) :: ixi^
l, ixo^
l
263 double precision,
intent(in) :: x(ixi^s,1:
ndim)
264 double precision,
intent(inout) :: w(ixi^s,1:nw)
265 end subroutine sts_before_first_cycle
267 subroutine sts_after_last_cycle(ixI^L, ixO^L, w, x)
269 integer,
intent(in) :: ixi^
l, ixo^
l
270 double precision,
intent(in) :: x(ixi^s,1:
ndim)
271 double precision,
intent(inout) :: w(ixi^s,1:nw)
272 end subroutine sts_after_last_cycle
285 subroutine sts_error_handling(w, x, ixI^L, ixO^L, step)
288 integer,
intent(in) :: ixi^
l,ixo^
l
289 double precision,
intent(inout) :: w(ixi^s,1:nw)
290 double precision,
intent(in) :: x(ixi^s,1:
ndim)
291 integer,
intent(in) :: step
292 end subroutine sts_error_handling
299 function sts_get_ncycles1(dt,dtnew,dt_modified)
result(is)
300 double precision,
intent(in) :: dtnew
301 double precision,
intent(inout) :: dt
302 logical,
intent(inout) :: dt_modified
305 double precision :: ss
308 ss = dtnew*((2.d0*sts_ncycles+1)**2-9.d0)/16.d0
316 if(ss .le. 1.d0)
then
319 is=ceiling((dsqrt(9.d0+16.d0*ss)-1.d0)*0.5d0)
324 end function sts_get_ncycles1
327 function sts_get_ncycles2(dt,dtnew,dt_modified)
result(is)
328 double precision,
intent(in) :: dtnew
329 double precision,
intent(inout) :: dt
330 logical,
intent(inout) :: dt_modified
333 double precision :: ss,rr
338 ncycles = sts_ncycles
340 ss=sum_chev(nu_sts,ncycles,rr)
349 end function sts_get_ncycles2
355 double precision,
intent(inout) :: my_dt
356 double precision :: my_dt1
357 logical :: dt_modified, dt_modified1, dt_modified2
359 double precision :: dtnew,dtmin_mype
360 double precision ::
dx^
d, ss
361 integer:: iigrid, igrid, ncycles
362 type(sts_term),
pointer :: temp,oldtemp
365 dt_modified = .false.
366 do while(
associated(temp))
367 dt_modified2 = .false.
370 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
377 dtmin_mype=min(dtmin_mype, sts_dtpar * temp%sts_getdt(ps(igrid)%w,ixg^
ll,
ixm^
ll,
dx^
d,ps(igrid)%x))
380 call mpi_allreduce(dtmin_mype,dtnew,1,mpi_double_precision,mpi_min,
icomm,
ierrmpi)
381 temp%s = sts_get_ncycles(my_dt,dtnew,dt_modified2)
389 if(dt_modified2)
then
394 dt_modified1 = .false.
395 do while(.not.
associated(oldtemp,temp))
396 oldtemp%s = sts_get_ncycles(my_dt1,oldtemp%dt_expl,dt_modified1)
398 if(dt_modified1)
call mpistop(
"sts dt modified twice")
399 oldtemp=>oldtemp%next
408 pure FUNCTION chev(j,nu,N)
411 double precision,
INTENT(IN) :: nu
412 INTEGER,
INTENT(IN) :: j, n
413 double precision :: chev
415 chev = 1d0 / ((-1d0 + nu)*cos(((2d0*j - 1d0) / n)* (
dpi/2d0)) + 1d0 + nu)
419 FUNCTION sum_chev(nu,N,limMax)
420 double precision,
intent(in) :: nu,limmax
421 integer,
intent(inout) :: n
422 double precision :: sum_chev, tmp
428 do while (j < n .and. sum_chev < limmax)
429 sum_chev = sum_chev + chev(j,nu,n)
433 END FUNCTION sum_chev
437 subroutine sts_add_source2(my_dt)
444 double precision,
intent(in) :: my_dt
445 double precision,
allocatable :: bj(:)
446 double precision :: sumbj,dtj
448 integer:: iigrid, igrid, j, ixc^
l
449 logical :: stagger_flag=.false., prolong_flag=.false., coarsen_flag=.false.
450 type(sts_term),
pointer :: temp
458 do while(
associated(temp))
460 if(.not.temp%evolve_magnetic_field)
then
471 if(
associated(temp%sts_before_first_cycle))
then
476 do iigrid=1,igridstail; igrid=igrids(iigrid);
479 call temp%sts_before_first_cycle(ixg^
ll,ixg^
ll,ps(igrid)%w,ps(igrid)%x)
483 allocate(bj(1:temp%s))
485 bj(j) = chev(j,nu_sts,sts_ncycles)
495 if(.not. temp%types_initialized)
then
497 if(temp%nflux>temp%nwbc)
then
513 temp%types_initialized = .true.
518 if(j .eq. temp%s .and. (sumbj + bj(j)) * temp%dt_expl > my_dt)
then
519 dtj = my_dt - sumbj * temp%dt_expl
521 dtj = bj(j)* temp%dt_expl
523 sumbj = sumbj + bj(j)
526 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
529 call temp%sts_set_sources(ixg^
ll,
ixm^
ll,ps(igrid)%w,ps(igrid)%x,ps1(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
530 if(temp%nflux>
ndir)
then
531 ps(igrid)%w(
ixm^t,temp%startVar)=ps(igrid)%w(
ixm^t,temp%startVar)+dtj*ps1(igrid)%w(
ixm^t,temp%startVar)
533 ps(igrid)%ws(ixc^s,1:nws)=ps(igrid)%ws(ixc^s,1:nws)+dtj*ps1(igrid)%w(ixc^s,iw_mag(1:nws))
539 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
542 call temp%sts_set_sources(ixg^
ll,
ixm^
ll,ps(igrid)%w,ps(igrid)%x,ps1(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
543 ps(igrid)%w(
ixm^t,temp%startVar:temp%endVar)=ps(igrid)%w(
ixm^t,temp%startVar:temp%endVar)+&
544 dtj*ps1(igrid)%w(
ixm^t,temp%startVar:temp%endVar)
549 if(fix_conserve_at_step)
then
557 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
563 if(
associated(temp%sts_handle_errors))
then
565 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
568 call temp%sts_handle_errors(ps(igrid)%w,ps(igrid)%x,ixg^
ll,
ixm^
ll,j)
573 if(temp%nflux>temp%nwbc.and.temp%s==j)
then
587 if(
associated(temp%sts_after_last_cycle))
then
588 do iigrid=1,igridstail; igrid=igrids(iigrid);
591 call temp%sts_after_last_cycle(ixg^
ll,ixg^
ll,ps(igrid)%w,ps(igrid)%x)
598 if(.not.temp%evolve_magnetic_field)
then
621 do iigrid=1,igridstail; igrid=igrids(iigrid);
627 end subroutine sts_add_source2
631 subroutine sts_add_source1(my_dt)
639 double precision,
intent(in) :: my_dt
640 double precision :: dtj
641 double precision :: omega1,cmu,cmut,cnu,cnut,one_mu_nu
642 double precision,
allocatable :: bj(:)
643 integer:: iigrid, igrid, j, ixc^
l, ixgext^
l
644 logical :: evenstep, stagger_flag=.false., prolong_flag=.false., coarsen_flag=.false., total_energy_flag=.true.
645 type(sts_term),
pointer :: temp
646 type(state),
dimension(:),
pointer :: tmpps1, tmpps2
654 do while(
associated(temp))
656 if(.not.temp%evolve_magnetic_field)
then
667 if(
associated(temp%sts_before_first_cycle))
then
675 do iigrid=1,igridstail; igrid=igrids(iigrid);
678 call temp%sts_before_first_cycle(ixg^
ll,ixg^
ll,ps(igrid)%w,ps(igrid)%x)
679 if(.not.
allocated(ps2(igrid)%w))
allocate(ps2(igrid)%w(ixg^t,1:nw))
680 if(.not.
allocated(ps3(igrid)%w))
allocate(ps3(igrid)%w(ixg^t,1:nw))
681 if(.not.
allocated(ps4(igrid)%w))
allocate(ps4(igrid)%w(ixg^t,1:nw))
682 ps1(igrid)%w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
683 ps2(igrid)%w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
688 ixgext^
l=ixg^
ll^ladd1;
690 do iigrid=1,igridstail; igrid=igrids(iigrid);
691 if(.not.
allocated(ps2(igrid)%w))
then
694 if(.not.
allocated(ps3(igrid)%w))
allocate(ps3(igrid)%w(ixg^t,1:nw))
695 if(.not.
allocated(ps4(igrid)%w))
allocate(ps4(igrid)%w(ixg^t,1:nw))
696 ps1(igrid)%w=ps(igrid)%w
697 ps2(igrid)%w=ps(igrid)%w
698 ps1(igrid)%ws=ps(igrid)%ws
699 ps2(igrid)%ws=ps(igrid)%ws
704 do iigrid=1,igridstail; igrid=igrids(iigrid);
705 if(.not.
allocated(ps2(igrid)%w))
allocate(ps2(igrid)%w(ixg^t,1:nw))
706 if(.not.
allocated(ps3(igrid)%w))
allocate(ps3(igrid)%w(ixg^t,1:nw))
707 if(.not.
allocated(ps4(igrid)%w))
allocate(ps4(igrid)%w(ixg^t,1:nw))
708 ps1(igrid)%w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
709 ps2(igrid)%w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
715 allocate(bj(0:temp%s))
719 omega1=4.d0/dble(temp%s**2+temp%s-2)
733 if(.not. temp%types_initialized)
then
735 if(temp%nflux>temp%nwbc)
then
751 temp%types_initialized = .true.
756 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
760 call temp%sts_set_sources(ixg^
ll,
ixm^
ll,ps(igrid)%w,ps(igrid)%x,ps4(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
762 ps3(igrid)%w(ixc^s,temp%startVar:temp%endVar) = my_dt * ps4(igrid)%w(ixc^s,temp%startVar:temp%endVar)
763 if(temp%nflux>
ndir)
then
764 ps1(igrid)%w(
ixm^t,temp%startVar) = ps1(igrid)%w(
ixm^t,temp%startVar) + cmut * ps3(igrid)%w(
ixm^t,temp%startVar)
766 ps1(igrid)%ws(ixc^s,1:nws) = ps1(igrid)%ws(ixc^s,1:nws) + cmut * ps3(igrid)%w(ixc^s,iw_mag(1:nws))
772 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
775 call temp%sts_set_sources(ixg^
ll,
ixm^
ll,ps(igrid)%w,ps(igrid)%x,ps4(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
777 ps3(igrid)%w(
ixm^t,temp%startVar:temp%endVar) = my_dt * ps4(igrid)%w(
ixm^t,temp%startVar:temp%endVar)
778 ps1(igrid)%w(
ixm^t,temp%startVar:temp%endVar) = ps1(igrid)%w(
ixm^t,temp%startVar:temp%endVar) + &
779 cmut * ps3(igrid)%w(
ixm^t,temp%startVar:temp%endVar)
783 if(fix_conserve_at_step)
then
791 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
798 if(
associated(temp%sts_handle_errors))
then
800 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
803 call temp%sts_handle_errors(ps1(igrid)%w,ps1(igrid)%x,ixg^
ll,
ixm^
ll,1)
807 if(temp%nflux>temp%nwbc.and.temp%s==1)
then
826 bj(j)=dble(j**2+j-2)/dble(2*j*(j+1))
827 cmu=dble(2*j-1)/dble(j)*bj(j)/bj(j-1)
829 cnu=dble(1-j)/dble(j)*bj(j)/bj(j-2)
830 cnut=(bj(j-1)-1.d0)*cmut
831 one_mu_nu=1.d0-cmu-cnu
843 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
849 call temp%sts_set_sources(ixg^
ll,
ixm^
ll,tmpps1(igrid)%w,ps(igrid)%x,ps4(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
850 if(temp%nflux>
ndir)
then
851 tmpps2(igrid)%w(
ixm^t,temp%startVar)=cmu*tmpps1(igrid)%w(
ixm^t,temp%startVar)+&
852 cnu*tmpps2(igrid)%w(
ixm^t,temp%startVar)+one_mu_nu*ps(igrid)%w(
ixm^t,temp%startVar)+&
853 dtj*ps4(igrid)%w(
ixm^t,temp%startVar)+cnut*ps3(igrid)%w(
ixm^t,temp%startVar)
855 tmpps2(igrid)%ws(ixc^s,1:nws)=cmu*tmpps1(igrid)%ws(ixc^s,1:nws)+&
856 cnu*tmpps2(igrid)%ws(ixc^s,1:nws)+one_mu_nu*ps(igrid)%ws(ixc^s,1:nws)+&
857 dtj*ps4(igrid)%w(ixc^s,iw_mag(1:nws))+cnut*ps3(igrid)%w(ixc^s,iw_mag(1:nws))
863 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
869 call temp%sts_set_sources(ixg^
ll,
ixm^
ll,tmpps1(igrid)%w,ps(igrid)%x,ps4(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
870 tmpps2(igrid)%w(
ixm^t,temp%startVar:temp%endVar)=cmu*tmpps1(igrid)%w(
ixm^t,temp%startVar:temp%endVar)+&
871 cnu*tmpps2(igrid)%w(
ixm^t,temp%startVar:temp%endVar)+one_mu_nu*ps(igrid)%w(
ixm^t,temp%startVar:temp%endVar)+&
872 dtj*ps4(igrid)%w(
ixm^t,temp%startVar:temp%endVar)+cnut*ps3(igrid)%w(
ixm^t,temp%startVar:temp%endVar)
876 if(fix_conserve_at_step)
then
884 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
890 if(
associated(temp%sts_handle_errors))
then
892 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
895 call temp%sts_handle_errors(tmpps2(igrid)%w,ps(igrid)%x,ixg^
ll,
ixm^
ll,j)
899 if(temp%nflux>temp%nwbc.and.temp%s==j)
then
911 evenstep=.not.evenstep
914 if(
associated(temp%sts_after_last_cycle))
then
916 do iigrid=1,igridstail; igrid=igrids(iigrid);
919 ps(igrid)%w(ixg^t,temp%startVar:temp%endVar)=tmpps2(igrid)%w(ixg^t,temp%startVar:temp%endVar)
920 call temp%sts_after_last_cycle(ixg^
ll,ixg^
ll,ps(igrid)%w,ps(igrid)%x)
929 do iigrid=1,igridstail; igrid=igrids(iigrid);
930 ps(igrid)%w(ixg^t,temp%startVar:temp%endVar)=tmpps2(igrid)%w(ixg^t,temp%startVar:temp%endVar)
931 ps(igrid)%ws=tmpps2(igrid)%ws
936 do iigrid=1,igridstail; igrid=igrids(iigrid);
937 ps(igrid)%w(ixg^t,temp%startVar:temp%endVar)=tmpps2(igrid)%w(ixg^t,temp%startVar:temp%endVar)
945 if(.not.temp%evolve_magnetic_field)
then
968 do iigrid=1,igridstail; igrid=igrids(iigrid);
974 end subroutine sts_add_source1
subroutine, public alloc_state(igrid, s, ixgl, ixgextl, alloc_once_for_ps)
allocate memory to physical state of igrid node
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter dpi
Pi.
Module for flux conservation near refinement boundaries.
subroutine, public init_comm_fix_conserve(idimlim, nwfluxin)
subroutine, public fix_edges(psuse, idimlim)
subroutine, public recvflux(idimlim)
subroutine, public sendflux(idimlim)
subroutine, public fix_conserve(psb, idimlim, nw0, nwfluxin)
Module with geometry-related routines (e.g., divergence, curl)
update ghost cells of all blocks including physical boundaries
integer, dimension( :^d &), pointer type_recv_r
integer, dimension(^nd, 0:3) l
subroutine getbc(time, qdt, psb, nwstart, nwbc)
do update ghost cells of all blocks including physical boundaries
integer, dimension(0:3^d &), target type_send_p_f
integer, dimension( :^d &), pointer type_send_p
integer, dimension( :^d &), pointer type_send_srl
subroutine create_bc_mpi_datatype(nwstart, nwbc)
integer, dimension(0:3^d &), target type_recv_r_f
integer, dimension(-1:1^d &), target type_recv_srl_f
integer, dimension(-1:1^d &), target type_send_r_f
integer, dimension(-1:1^d &), target type_send_srl_f
integer, dimension( :^d &), pointer type_send_r
integer, dimension( :^d &), pointer type_recv_p
integer, dimension(0:3^d &), target type_recv_p_f
integer, dimension( :^d &), pointer type_recv_srl
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
integer, parameter unitpar
file handle for IO
double precision global_time
The global simulation time.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer icomm
The MPI communicator.
logical coarsenprimitive
coarsen primitive variables in level-jump ghost cells
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
logical time_advance
do time evolving
logical prolongprimitive
prolongate primitive variables in level-jump ghost cells
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
This module defines the procedures of a physics module. It contains function pointers for the various...
logical phys_partial_ionization
Solve partially ionized one-fluid plasma.
logical phys_total_energy
Solve total energy equation or not.
procedure(sub_update_temperature), pointer phys_update_temperature
procedure(sub_face_to_center), pointer phys_face_to_center
Module for handling problematic values in simulations, such as negative pressures.
Generic supertimestepping method which can be used for multiple source terms in the governing equatio...
integer, public sourcetype_sts
pure logical function, public is_sts_initialized()
integer, parameter, public sourcetype_sts_prior
logical function, public set_dt_sts_ncycles(my_dt)
This sets the explicit dt and calculates the number of cycles for each of the terms implemented with ...
subroutine, public add_sts_method(sts_getdt, sts_set_sources, startvar, nflux, startwbc, nwbc, evolve_b)
subroutine which added programatically a term to be calculated using STS Params: sts_getdt function c...
type(sts_term), pointer head_sts_terms
subroutine, public set_conversion_methods_to_head(sts_before_first_cycle, sts_after_last_cycle)
Set the hooks called before the first cycle and after the last cycle in the STS update This method sh...
integer, parameter, public sourcetype_sts_split
subroutine, public set_error_handling_to_head(sts_error_handling)
Set the hook of error handling in the STS update. This method is called before updating the BC....
subroutine, public sts_init()
Initialize sts module.
integer, parameter, public sourcetype_sts_after
procedure(subr3), pointer, public sts_add_source