MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_point_searching.t
Go to the documentation of this file.
1 ! for Cartesian coordinate system
4  use mod_physics
6  implicit none
7 
8 contains
9 
10  subroutine get_point_w_ingrid(igrid,xp,wp,variable_type)
11  double precision :: xp(1:ndim),wp(1:nw)
12  character(*) :: variable_type
13  integer, intent(in) :: igrid
14 
15  double precision :: dxb^D,xd^D,xb^L,temp
16  integer :: ixI^L,ixO^L
17  integer :: ixbl^D,ix^D,ixA^L,j,ingrid
18  double precision :: factor(0:1^D&)
19 
20  ixi^l=ixg^ll;
21  ixo^l=ixm^ll;
22  ^d&xbmin^d=rnode(rpxmin^d_,igrid)-rnode(rpdx^d_,igrid);
23  ^d&xbmax^d=rnode(rpxmax^d_,igrid)+rnode(rpdx^d_,igrid);
24  ingrid=0
25  {if (xp(^db)>xbmin^db .and. xp(^db)<xbmax^db) ingrid=ingrid+1\}
26 
27  if (ingrid==ndim) then
28  ^d&dxb^d=rnode(rpdx^d_,igrid)\
29  ^d&ixbl^d=floor((xp(^d)-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d\
30  ^d&xd^d=(xp(^d)-ps(igrid)%x(ixbl^dd,^d))/dxb^d\
31  ^d&ixamin^d=ixbl^d\
32  ^d&ixamax^d=ixbl^d+1\
33 
34  {do ix^d=0,1\}
35  factor(ix^d)={abs(1-ix^d-xd^d)*}
36  {enddo\}
37 
38  select case(variable_type)
39  case('primitive')
40  call phys_to_primitive(ixi^l,ixa^l,ps(igrid)%w,ps(igrid)%x)
41  case('conserved')
42 
43  case default
44  call mpistop("get_point_w_ingrid: Unknown variable type!")
45  end select
46 
47  wp=0.d0
48  {do ix^d=0,1\}
49  do j=1,nw
50  wp(j)=wp(j)+factor(ix^d)*ps(igrid)%w(ixbl^d+ix^d,j)
51  enddo
52  {enddo\}
53 
54  select case(variable_type)
55  case('primitive')
56  call phys_to_conserved(ixi^l,ixa^l,ps(igrid)%w,ps(igrid)%x)
57  case('conserved')
58 
59  case default
60  call mpistop("get_point_w_ingrid: Unknown variable type!")
61  end select
62 
63  if(physics_type=='mhd') then
64  wp(iw_mag(1):iw_mag(ndir))=0.d0
65  {do ix^d=0,1\}
66  do j=1,ndir
67  if (b0field) then
68  temp=ps(igrid)%w(ixbl^d+ix^d,iw_mag(j))+&
69  ps(igrid)%B0(ixbl^d+ix^d,j,0)
70  wp(iw_mag(j))=wp(iw_mag(j))+factor(ix^d)*temp
71  else
72  temp=ps(igrid)%w(ixbl^d+ix^d,iw_mag(j))
73  wp(iw_mag(j))=wp(iw_mag(j))+factor(ix^d)*temp
74  endif
75  enddo
76  {enddo\}
77  endif
78  else
79  call mpistop("get_point_w_ingrid: The point is not in given grid!")
80  endif
81 
82  end subroutine get_point_w_ingrid
83 
84  subroutine get_point_w(xp,wp,variable_type)
85  ! for given point (xp), provide the plasma parameters (wp) at this point
86 
87  double precision :: xp(1:ndim),wp(1:nw)
88  character(*) :: variable_type
89 
90  double precision :: x3d(3)
91  integer :: indomain,ipe,igrid,j
92 
93  indomain=0
94  {if (xp(^db)>=xprobmin^db .and. xp(^db)<xprobmax^db) indomain=indomain+1\}
95  if (indomain==ndim) then
96 
97  ! find pe and igrid
98  x3d=0.d0
99  do j=1,ndim
100  x3d(j)=xp(j)
101  enddo
102  call find_particle_ipe(x3d,igrid,ipe)
103 
104  !do interpolation to get point values
105  if (mype==ipe) then
106  call get_point_w_ingrid(igrid,xp,wp,variable_type)
107  endif
108 
109  call mpi_bcast(wp,nw,mpi_double_precision,ipe,icomm,ierrmpi)
110  endif
111 
112  end subroutine get_point_w
113 
114  subroutine get_cell_w(xp,wc,variable_type)
115  ! for given point (xp), looking for corresponding cell and then provide
116  ! the plasma parameters (wc) of this cell
117 
118  double precision :: xp(1:ndim),wc(1:nw)
119  character(*) :: variable_type
120 
121  double precision :: x3d(3)
122  double precision :: dxb^D,xb^L
123  integer :: indomain,ixO^L,ixb^D,ixA^L,ixI^L
124  integer :: ipe,igrid,j
125 
126  indomain=0
127  {if (xp(^db)>=xprobmin^db .and. xp(^db)<xprobmax^db) indomain=indomain+1\}
128  if (indomain==ndim) then
129 
130  ! find pe and igrid
131  x3d=0.d0
132  do j=1,ndim
133  x3d(j)=xp(j)
134  enddo
135  call find_particle_ipe(x3d,igrid,ipe)
136 
137  if (mype==ipe) then
138  ! looking for the cell index
139  ^d&xbmin^d=rnode(rpxmin^d_,igrid)\
140  ^d&xbmax^d=rnode(rpxmax^d_,igrid)\
141  ^d&dxb^d=rnode(rpdx^d_,igrid)\
142  ^d&ixomin^d=ixmlo^d\
143  ^d&iximin^d=ixglo^d\
144  ^d&iximax^d=ixghi^d\
145  ^d&ixb^d=floor((xp(^d)-xbmin^d)/dxb^d)+ixomin^d\
146  ^d&ixamin^d=ixb^d\
147  ^d&ixamax^d=ixb^d\
148 
149  wc=0.d0
150  select case(variable_type)
151  case('primitive')
152  call phys_to_primitive(ixi^l,ixa^l,ps(igrid)%w,ps(igrid)%x)
153  do j=1,nw
154  wc(j)=ps(igrid)%w(ixb^d,j)
155  enddo
156  call phys_to_conserved(ixi^l,ixa^l,ps(igrid)%w,ps(igrid)%x)
157  case('conserved')
158  do j=1,nw
159  wc(j)=ps(igrid)%w(ixb^d,j)
160  enddo
161  case default
162  call mpistop("get_point_w: Unknown variable type!")
163  end select
164 
165  ! for magnetic field
166  if(physics_type=='mhd') then
167  wc(iw_mag(1):iw_mag(ndir))=0.d0
168  do j=1,ndir
169  if (b0field) then
170  wc(iw_mag(j))=ps(igrid)%w(ixb^d,iw_mag(j))+&
171  ps(igrid)%B0(ixb^d,j,0)
172  else
173  wc(iw_mag(j))=ps(igrid)%w(ixb^d,iw_mag(j))
174  endif
175  enddo
176  endif
177  endif
178 
179  call mpi_bcast(wc,nw,mpi_double_precision,ipe,icomm,ierrmpi)
180  endif
181 
182  end subroutine get_cell_w
183 
184  subroutine get_cell_index(xp,ipe,igrid,ixc^D)
185  ! for given point (xp), provide the corrsponding processor (ipe),
186  ! grid number (igrid) and cell index (ixc^D)
187 
188  double precision :: xp(1:ndim)
189  integer :: ipe,igrid,ixc^D
190 
191  double precision :: x3d(1:3)
192  double precision :: dxb^D,xb^L
193  integer :: indomain,ixO^L,j
194  integer :: datas(1:ndim+2)
195 
196  indomain=0
197  {if (xp(^db)>=xprobmin^db .and. xp(^db)<xprobmax^db) indomain=indomain+1\}
198 
199  if (indomain==ndim) then
200  ! find pe and igrid
201  x3d=0.d0
202  do j=1,ndim
203  x3d(j)=xp(j)
204  enddo
205  call find_particle_ipe(x3d,igrid,ipe)
206 
207  if (mype==ipe) then
208  ! looking for the cell index
209  ^d&xbmin^d=rnode(rpxmin^d_,igrid)\
210  ^d&xbmax^d=rnode(rpxmax^d_,igrid)\
211  ^d&dxb^d=rnode(rpdx^d_,igrid)\
212  ^d&ixomin^d=ixmlo^d\
213  ^d&ixc^d=floor((xp(^d)-xbmin^d)/dxb^d)+ixomin^d\
214  ^d&datas(^d)=ixc^d\
215  datas(ndim+1)=ipe
216  datas(ndim+2)=igrid
217  endif
218 
219  call mpi_bcast(datas,ndim+2,mpi_integer,ipe,icomm,ierrmpi)
220  endif
221 
222  ^d&ixc^d=datas(^d)\
223  ipe=datas(ndim+1)
224  igrid=datas(ndim+2)
225 
226  end subroutine get_cell_index
227 
228 end module mod_point_searching
subroutine mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: comm_lib.t:194
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer ixghi
Upper index of grid block arrays.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer ndir
Number of spatial dimensions (components) for vector variables.
integer ixm
the mesh range (within a block with ghost cells)
integer ierrmpi
A global MPI error return code.
logical b0field
split magnetic field as background B0 field
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, parameter ixglo
Lower index of grid block arrays (always 1)
Module with shared functionality for all the particle movers.
subroutine find_particle_ipe(x, igrid_particle, ipe_particle)
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
procedure(sub_convert), pointer phys_to_conserved
Definition: mod_physics.t:58
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition: mod_physics.t:19
subroutine get_point_w(xp, wp, variable_type)
subroutine get_cell_index(xp, ipe, igrid, ixcD)
subroutine get_cell_w(xp, wc, variable_type)
subroutine get_point_w_ingrid(igrid, xp, wp, variable_type)