MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_source.t
Go to the documentation of this file.
1 !> Module for handling split source terms (split from the fluxes)
2 module mod_source
3  implicit none
4  public
5 
6  !> How to apply dimensional splitting to the source terms, see
7  !> @ref discretization.md
8  integer :: sourcesplit =-1
9  integer, parameter :: sourcesplit_sfs = 0
10  integer, parameter :: sourcesplit_sf = 1
11  integer, parameter :: sourcesplit_ssf = 2
12  integer, parameter :: sourcesplit_ssfss = 3
13 
14  public :: add_split_source
15  public :: addsource2
16 
17 contains
18 
19  subroutine add_split_source(prior)
25  use mod_comm_lib, only: mpistop
26 
27  logical, intent(in) :: prior
28  ! This variable, later allocated on the thread stack, causes segmentation fault
29  ! when openmp is used with intel. That could be solved otherwise, by increasing
30  ! the thread stack size, but not using it at all could speed up.
31  double precision :: w1(ixg^t,1:nw)
32  double precision :: qt
33  integer :: iigrid, igrid
34  logical :: src_active
35 
36  ! add stiff source terms via super time stepping
37  if(is_sts_initialized()) then
38  select case (sourcetype_sts)
40  if(prior) then
41  call sts_add_source(dt)
42  end if
44  if(.not. prior) then
45  call sts_add_source(dt)
46  end if
48  call sts_add_source(0.5d0*dt)
49  end select
50  end if
51  src_active = .false.
52 
53  if ((.not.prior).and.&
55 
56  if (prior) then
57  qt=global_time
58  else
59  qt=global_time+dt
60  end if
61 
62  if(any_source_split) then
63  ! add normal split source terms
64  select case (sourcesplit)
65  case (sourcesplit_sfs)
66  !$OMP PARALLEL DO PRIVATE(igrid)
67  do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
68  block=>ps(igrid)
69  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
70  w1=ps(igrid)%w
71  call phys_to_primitive(ixg^ll,ixg^ll,w1,ps(igrid)%x)
72  call addsource2(0.5d0*dt,0.5d0,ixg^ll,ixm^ll,1,nw,qt,ps(igrid)%w,w1,qt,ps(igrid)%w,&
73  ps(igrid)%x,.true.,src_active)
74  end do
75  !$OMP END PARALLEL DO
76  case (sourcesplit_sf)
77  !$OMP PARALLEL DO PRIVATE(igrid)
78  do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
79  block=>ps(igrid)
80  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
81  w1=ps(igrid)%w
82  call phys_to_primitive(ixg^ll,ixg^ll,w1,ps(igrid)%x)
83  call addsource2(dt ,1d0,ixg^ll,ixm^ll,1,nw,qt,ps(igrid)%w,w1,qt,ps(igrid)%w,&
84  ps(igrid)%x,.true.,src_active)
85  end do
86  !$OMP END PARALLEL DO
87  case (sourcesplit_ssf)
88  !$OMP PARALLEL DO PRIVATE(igrid)
89  do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
90  block=>ps(igrid)
91  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
92  w1=ps(igrid)%w
93  call phys_to_primitive(ixg^ll,ixg^ll,w1,ps(igrid)%x)
94  call addsource2(0.5d0*dt,0.5d0,ixg^ll,ixg^ll,1,nw,qt,ps(igrid)%w,w1,qt,ps(igrid)%w,&
95  ps(igrid)%x,.true.,src_active)
96  call addsource2(dt ,1d0,ixg^ll,ixm^ll,1,nw,qt,ps(igrid)%w,w1,qt,ps(igrid)%w,&
97  ps(igrid)%x,.true.,src_active)
98  end do
99  !$OMP END PARALLEL DO
100  case (sourcesplit_ssfss)
101  !$OMP PARALLEL DO PRIVATE(igrid)
102  do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
103  block=>ps(igrid)
104  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
105  w1=ps(igrid)%w
106  call phys_to_primitive(ixg^ll,ixg^ll,w1,ps(igrid)%x)
107  call addsource2(0.25d0*dt,0.25d0,ixg^ll,ixg^ll,1,nw,qt,ps(igrid)%w,w1,qt,ps(igrid)%w,&
108  ps(igrid)%x,.true.,src_active)
109  call addsource2(0.5d0*dt,0.5d0,ixg^ll,ixm^ll,1,nw,qt,ps(igrid)%w,w1,qt,ps(igrid)%w,&
110  ps(igrid)%x,.true.,src_active)
111  end do
112  !$OMP END PARALLEL DO
113  case default
114  write(unitterm,*)'No such type of sourcesplit=',sourcesplit
115  call mpistop("Error: Unknown type of sourcesplit!")
116  end select
117  end if
118 
119  if (.not. prior .and. associated(phys_global_source_after)) then
120  call phys_global_source_after(dt, qt, src_active)
121  end if
122 
123  if (src_active) then
124  call getbc(qt,0.d0,ps,iwstart,nwgc,phys_req_diagonal)
125  end if
126 
127  end subroutine add_split_source
128 
129  !> Add source within ixO for iws: w=w+qdt*S[wCT]
130  subroutine addsource2(qdt,dtfactor,ixI^L,ixO^L,iw^LIM,qtC,wCT,wCTprim,qt,&
131  w,x,qsourcesplit,src_active)
133  use mod_physics, only: phys_add_source
134  use mod_usr_methods, only: usr_source
135  ! differences with VAC is in iw^LIM and in declaration of ranges for wCT,w
136 
137  integer, intent(in) :: ixi^l, ixo^l, iw^lim
138  double precision, intent(in) :: qdt, dtfactor, qtc, qt
139  double precision, intent(in) :: wct(ixi^s,1:nw), wctprim(ixi^s,1:nw), x(ixi^s,1:ndim)
140  double precision, intent(inout) :: w(ixi^s,1:nw)
141  logical, intent(in) :: qsourcesplit
142  logical, intent(inout), optional :: src_active
143 
144  logical :: tmp_active
145 
146  tmp_active = .false.
147 
148  ! user defined sources, typically explicitly added
149  if ((qsourcesplit .eqv. source_split_usr) .and. associated(usr_source)) then
150  tmp_active = .true.
151  call usr_source(qdt,ixi^l,ixo^l,iw^lim,qtc,wct,qt,w,x)
152  end if
153 
154  ! physics defined sources, typically explicitly added,
155  ! along with geometrical source additions
156  call phys_add_source(qdt,dtfactor,ixi^l,ixo^l,wct,wctprim,w,x,qsourcesplit,tmp_active)
157 
158  if (present(src_active)) src_active = src_active .or. tmp_active
159 
160  end subroutine addsource2
161 
162 end module mod_source
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
update ghost cells of all blocks including physical boundaries
subroutine getbc(time, qdt, psb, nwstart, nwbc, req_diag)
do update ghost cells of all blocks including physical boundaries
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.
logical source_split_usr
Use split or unsplit way to add user's source terms, default: unsplit.
logical any_source_split
if any normal source term is added in split fasion
double precision global_time
The global simulation time.
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer, dimension(:), allocatable, parameter d
integer, parameter unitterm
Unit for standard output.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision, dimension(ndim) dxlevel
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition: mod_physics.t:59
logical phys_req_diagonal
Whether the physics routines require diagonal ghost cells, for example for computing a curl.
Definition: mod_physics.t:24
procedure(sub_global_source), pointer phys_global_source_after
Definition: mod_physics.t:72
procedure(sub_add_source), pointer phys_add_source
Definition: mod_physics.t:71
Module for handling split source terms (split from the fluxes)
Definition: mod_source.t:2
subroutine, public addsource2(qdt, dtfactor, ixIL, ixOL, iwLIM, qtC, wCT, wCTprim, qt, w, x, qsourcesplit, src_active)
Add source within ixO for iws: w=w+qdt*S[wCT].
Definition: mod_source.t:132
integer sourcesplit
How to apply dimensional splitting to the source terms, see Discretization.
Definition: mod_source.t:8
integer, parameter sourcesplit_ssf
Definition: mod_source.t:11
integer, parameter sourcesplit_ssfss
Definition: mod_source.t:12
integer, parameter sourcesplit_sfs
Definition: mod_source.t:9
integer, parameter sourcesplit_sf
Definition: mod_source.t:10
subroutine, public add_split_source(prior)
Definition: mod_source.t:20
Generic supertimestepping method 1) in amrvac.par in sts_list set the following parameters which have...
integer, public sourcetype_sts
pure logical function, public is_sts_initialized()
integer, parameter, public sourcetype_sts_prior
integer, parameter, public sourcetype_sts_split
integer, parameter, public sourcetype_sts_after
procedure(subr3), pointer, public sts_add_source
Module with all the methods that users can customize in AMRVAC.
procedure(source), pointer usr_source