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