8 subroutine trace_field_multi(xfm,wPm,wLm,dL,numL,numP,nwP,nwL,forwardm,ftype,tcondi)
36 integer,
intent(in) :: numL,numP,nwP,nwL
37 double precision,
intent(inout) :: xfm(numL,numP,ndim),wPm(numL,numP,nwP),wLm(numL,1+nwL)
38 double precision,
intent(in) :: dL
39 logical,
intent(in) :: forwardm(numL)
40 character(len=std_len),
intent(in) :: ftype,tcondi
42 double precision :: x3d(3),statusF(4+ndim),statusL(numL,4+ndim+nwL),statusS(numL,4+ndim+nwL)
43 double precision :: xf(numP,ndim),wP(numP,nwP),wL(1+nwL)
44 double precision,
allocatable :: data_send(:,:,:),data_recv(:,:,:)
45 integer :: indomain,ipe_now,igrid_now,igrid,j,iL
46 integer :: ipoint_in,ipoint_out,numSend,nRT,nRTmax
47 integer :: ipointm(numL),igridm(numL)
48 logical :: continueL(numL),myL(numL)
49 logical :: stopT,forward
51 if (tcondi/=
'TRAC')
then
57 xfm(1:numl,2:nump,:)=zero
67 {
if (xfm(il,1,^db)>=xprobmin^db .and. xfm(il,1,^db)<xprobmax^db) indomain=indomain+1\}
68 if (indomain==ndim)
then
69 if (tcondi/=
'TRAC') wlm(il,1)=1
80 if (
mype==ipe_now)
then
91 do while(stopt .eqv. .false.)
95 if (myl(il) .and. continuel(il))
then
98 xf(ipoint_in,:)=xfm(il,ipoint_in,:)
102 call find_points_in_pe(igrid,ipoint_in,xf,wp,wl,dl,nump,nwp,nwl,forward,ftype,tcondi,statusf)
103 ipoint_out=int(statusf(1))
104 xfm(il,ipoint_in:ipoint_out-1,:)=xf(ipoint_in:ipoint_out-1,:)
105 wpm(il,ipoint_in:ipoint_out-1,:)=wp(ipoint_in:ipoint_out-1,:)
115 statuss(il,1:4+ndim)=statusf(1:4+ndim)
116 statuss(il,4+ndim+1:4+ndim+nwl)=wl(2:1+nwl)
117 if (tcondi==
'TRAC') wlm(il,1)=ipoint_out-1
122 numsend=numl*(4+ndim+nwl)
123 call mpi_allreduce(statuss,statusl,numsend,mpi_double_precision,&
130 if (continuel(il))
then
131 ipointm(il)=int(statusl(il,1))
132 if (
mype==int(statusl(il,2))) myl(il)=.true.
133 igridm(il)=int(statusl(il,3))
134 if (int(statusl(il,4))==0)
then
137 continuel(il)=.false.
139 if (myl(il)) xfm(il,ipointm(il),1:ndim)=statusl(il,4+1:4+ndim)
140 if (tcondi/=
'TRAC') wlm(il,1)=ipointm(il)-1
141 wlm(il,2:1+nwl)=statusl(il,4+ndim+1:4+ndim+nwl)
147 if (tcondi/=
'TRAC')
then
150 if (nrtmax<int(wlm(il,1))) nrtmax=int(wlm(il,1))
152 numsend=numl*nrtmax*(ndim+nwp)
154 allocate(data_send(numl,nrtmax,ndim+nwp),data_recv(numl,nrtmax,ndim+nwp))
155 data_send(:,:,:)=zero
158 data_send(il,1:nrt,1:ndim)=xfm(il,1:nrt,1:ndim)
159 if (nwp>0) data_send(il,1:nrt,1+ndim:ndim+nwp)=wpm(il,1:nrt,1:nwp)
161 call mpi_allreduce(data_send,data_recv,numsend,mpi_double_precision,&
165 xfm(il,1:nrt,1:ndim)=data_recv(il,1:nrt,1:ndim)
166 if (nwp>0) wpm(il,1:nrt,1:nwp)=data_recv(il,1:nrt,1+ndim:ndim+nwp)
168 deallocate(data_send,data_recv)
201 integer,
intent(in) :: numP,nwP,nwL
202 double precision,
intent(inout) :: xf(numP,ndim),wP(numP,nwP),wL(1+nwL)
203 double precision,
intent(in) :: dL
204 logical,
intent(in) :: forward
205 character(len=std_len),
intent(in) :: ftype,tcondi
207 double precision :: x3d(3),statusF(4+ndim),status_bcast(4+ndim+nwL)
208 double precision,
allocatable :: data_send(:,:),data_recv(:,:)
209 integer :: indomain,ipoint_in,ipe_now,igrid_now,igrid,j
210 integer :: ipoint_out,ipe_next,igrid_next,numRT
221 {
if (xf(1,^db)>=xprobmin^db .and. xf(1,^db)<xprobmax^db) indomain=indomain+1\}
222 if (indomain==ndim)
then
233 if (
mype/=ipe_now) xf(1,:)=zero
236 call mpistop(
'Field tracing error: given point is not in simulation box!')
242 do while(stopt .eqv. .false.)
244 if (
mype==ipe_now)
then
247 call find_points_in_pe(igrid,ipoint_in,xf,wp,wl,dl,nump,nwp,nwl,forward,ftype,tcondi,statusf)
248 status_bcast(1:4+ndim)=statusf(1:4+ndim)
249 status_bcast(4+ndim+1:4+ndim+nwl)=wl(2:1+nwl)
252 call mpi_bcast(status_bcast,4+ndim+nwl,mpi_double_precision,ipe_now,
icomm,
ierrmpi)
253 statusf(1:4+ndim)=status_bcast(1:4+ndim)
254 wl(2:1+nwl)=status_bcast(4+ndim+1:4+ndim+nwl)
257 ipoint_out=int(statusf(1))
258 ipe_next=int(statusf(2))
259 igrid_next=int(statusf(3))
260 if (int(statusf(4))==1)
then
264 if (
mype==ipe_next)
then
266 xf(ipoint_out,j)=statusf(4+j)
269 xf(ipoint_out,:)=zero
278 if (tcondi/=
'TRAC')
then
280 allocate(data_send(numrt,ndim+nwp),data_recv(numrt,ndim+nwp))
283 data_send(1:numrt,1:ndim)=xf(1:numrt,1:ndim)
284 if (nwp>0) data_send(1:numrt,1+ndim:ndim+nwp)=wp(1:numrt,1:nwp)
285 call mpi_allreduce(data_send,data_recv,numrt*(ndim+nwp),mpi_double_precision,&
287 xf(1:numrt,1:ndim)=data_recv(1:numrt,1:ndim)
288 if (nwp>0) wp(1:numrt,1:nwp)=data_recv(1:numrt,1+ndim:ndim+nwp)
289 deallocate(data_send,data_recv)
294 subroutine find_points_in_pe(igrid,ipoint_in,xf,wP,wL,dL,numP,nwP,nwL,forward,ftype,tcondi,statusF)
296 integer,
intent(inout) :: igrid
297 integer,
intent(in) :: ipoint_in,numP,nwP,nwL
298 double precision,
intent(inout) :: xf(numP,ndim),wP(numP,nwP),wL(1+nwL)
299 double precision,
intent(in) :: dL
300 logical,
intent(in) :: forward
301 character(len=std_len),
intent(in) :: ftype,tcondi
302 double precision,
intent(inout) :: statusF(4+ndim)
304 double precision :: xfout(ndim)
305 integer :: ipe_next,igrid_next,ip_in,ip_out,j,indomain
306 logical :: newpe,stopT
311 do while(newpe .eqv. .false.)
313 call find_points_interp(igrid,ip_in,ip_out,xf,wp,wl,nump,nwp,nwl,dl,forward,ftype,tcondi)
318 {
if (xf(ip_out,^db)>=xprobmin^db .and. xf(ip_out,^db)<=xprobmax^db) indomain=indomain+1\}
319 if (ip_out<nump .and. indomain==ndim)
then
320 if (tcondi/=
'TRAC')
then
342 statusf(3)=igrid_next
344 if (stopt) statusf(4)=1
346 statusf(4+j)=xf(ip_out,j)
350 if (newpe .eqv. .false.) igrid=igrid_next
361 integer,
intent(inout) :: igrid,igrid_next,ipe_next
362 double precision,
intent(in) :: xf1(ndim)
363 logical,
intent(inout) :: newpe,stopT
365 double precision :: dxb^D,xb^L,xbmid^D
366 double precision :: xbn^L
367 integer :: idn^D,my_neighbor_type,inblock
368 integer :: ic^D,inc^D,ipe_neighbor,igrid_neighbor
376 {
if (xf1(^d)<=xbmin^d) idn^d=-1\}
377 {
if (xf1(^d)>=xbmax^d) idn^d=1\}
378 my_neighbor_type=neighbor_type(idn^d,igrid)
379 igrid_neighbor=neighbor(1,idn^d,igrid)
380 ipe_neighbor=neighbor(2,idn^d,igrid)
383 select case(my_neighbor_type)
384 case (neighbor_boundary)
389 case(neighbor_coarse)
391 igrid_next=igrid_neighbor
392 ipe_next=ipe_neighbor
393 if (
mype==ipe_neighbor)
then
399 case(neighbor_sibling)
401 igrid_next=igrid_neighbor
402 ipe_next=ipe_neighbor
403 if (
mype==ipe_neighbor)
then
411 {xbmid^d=(xbmin^d+xbmax^d)/2.d0\}
413 {
if (xf1(^d)<=xbmin^d) inc^d=0\}
414 {
if (xf1(^d)>xbmin^d .and. xf1(^d)<=xbmid^d) inc^d=1\}
415 {
if (xf1(^d)>xbmid^d .and. xf1(^d)<xbmax^d) inc^d=2\}
416 {
if (xf1(^d)>=xbmax^d) inc^d=3\}
417 ipe_next=neighbor_child(2,inc^d,igrid)
418 igrid_next=neighbor_child(1,inc^d,igrid)
419 if (
mype==ipe_next)
then
428 subroutine find_points_interp(igrid,ip_in,ip_out,xf,wP,wL,numP,nwP,nwL,dL,forward,ftype,tcondi)
432 integer,
intent(in) :: igrid,ip_in,numP,nwP,nwL
433 integer,
intent(inout) :: ip_out
434 double precision,
intent(inout) :: xf(numP,ndim),wP(numP,nwP),wL(1+nwL)
435 double precision,
intent(in) :: dL
436 logical,
intent(in) :: forward
437 character(len=std_len),
intent(in) :: ftype,tcondi
439 double precision :: dxb^D,xb^L
440 double precision :: field(ixg^T,ndir)
441 double precision :: xs1(ndim),xs2(ndim),K1(ndim),K2(ndim)
442 double precision :: xfpre(ndim),xfnow(ndim),xfnext(ndim)
443 double precision :: Tpre,Tnow,Tnext,dTds,Lt,Lr,ds,T_bott,trac_delta
444 integer :: ip,inblock,ixI^L,ixO^L,j
452 if (tcondi/=
'TRAC')
then
461 mainloop:
do ip=ip_in,nump-1
465 call get_k(xs1,igrid,k1,ixi^l,dxb^d,ftype)
467 xs2(:)=xf(ip,:)+ds*k1(:)
469 xs2(:)=xf(ip,:)-ds*k1(:)
471 call get_k(xs2,igrid,k2,ixi^l,dxb^d,ftype)
473 xf(ip+1,:)=xf(ip,:)+ds*(0.5*k1(:)+0.5*k2(:))
475 xf(ip+1,:)=xf(ip,:)-ds*(0.5*k1(:)+0.5*k2(:))
480 if (tcondi/=
'TRAC')
then
482 call usr_set_field_w(igrid,ip,xf,wp,wl,nump,nwp,nwl,dl,forward,ftype,tcondi)
489 xfpre(:)=xf(ip,:)-ds*k1(:)
491 xfpre(:)=xf(ip,:)+ds*k1(:)
506 dtds=abs(tnext-tpre)/(2*ds)
517 if(lr>trac_delta*lt)
then
518 if (tnow>wl(2)) wl(2)=tnow
521 if (tnow>wl(3)) wl(3)=tnow
527 {
if (xf(ip+1,^db)>=xbmin^db .and. xf(ip+1,^db)<xbmax^db) inblock=inblock+1\}
529 if (inblock/=ndim)
exit mainloop
535 subroutine get_k(xfn,igrid,K,ixI^L,dxb^D,ftype)
538 integer :: ixI^L,igrid
539 double precision :: dxb^D
540 double precision :: xfn(ndim),K(ndim)
541 character(len=std_len) :: ftype
543 double precision :: xd^D
544 double precision :: field(0:1^D&,ndir),Fx(ndim),factor(0:1^D&)
545 double precision :: vector(ixI^S,1:ndir)
546 double precision :: Ftotal
547 integer :: ixb^D,ix^D,ixbl^D,ixO^L,j
549 ^d&ixbl^d=floor((xfn(^d)-ps(igrid)%x(iximin^dd,^d))/dxb^d)+iximin^d;
550 ^d&xd^d=(xfn(^d)-ps(igrid)%x(ixbl^dd,^d))/dxb^d;
552 ^d&ixomax^d=ixbl^d+1;
556 if(ftype==
'Bfield')
then
558 if(
allocated(iw_mag))
then
559 vector(ixo^s,1:ndir)=ps(igrid)%w(ixo^s,iw_mag(1:ndir))+ps(igrid)%B0(ixo^s,1:ndir,0)
561 vector(ixo^s,1:ndir)=ps(igrid)%B0(ixo^s,1:ndir,0)
564 vector(ixo^s,1:ndir)=ps(igrid)%w(ixo^s,iw_mag(1:ndir))
566 else if (ftype==
'Vfield')
then
567 call phys_get_v(ps(igrid)%w,ps(igrid)%x,ixi^l,ixo^l,vector)
570 factor(ix^d)={abs(1-ix^d-xd^d)*}
571 field(ix^d,1:ndir)=vector(ixbl^d+ix^d,1:ndir)
577 fx(j)=fx(j)+field(ix^d,j)*factor(ix^d)
581 if (ftype==
'Bfield' .or. ftype==
'Vfield')
then
585 fx(j)=fx(j)+field(ix^d,j)*factor(ix^d)
591 call mpistop(
'Field tracing error: wrong field type!')
596 ftotal=ftotal+(fx(j))**2
600 if (ftotal==0.d0)
then
604 k(1:ndim)=fx(1:ndim)/ftotal
612 integer,
intent(in) :: igrid,ixI^L
613 double precision,
intent(inout) :: xloc(ndim)
614 double precision,
intent(inout) :: Tloc
615 double precision,
intent(in) :: dxb^D
617 double precision :: xd^D
618 double precision :: factor(0:1^D&),Tnear(0:1^D&)
619 integer :: ixb^D,ix^D,ixbl^D,j,ixO^L
621 ^d&ixbl^d=floor((xloc(^d)-ps(igrid)%x(iximin^dd,^d))/dxb^d)+iximin^d;
622 ^d&xd^d=(xloc(^d)-ps(igrid)%x(ixbl^dd,^d))/dxb^d;
624 ^d&ixomax^d=ixomin^d+1;
627 factor(ix^d)={abs(1-ix^d-xd^d)*}
628 tnear(ix^d)=ps(igrid)%wextra(ixbl^d+ix^d,iw_tcoff)
634 tloc=tloc+tnear(ix^d)*factor(ix^d)
Module with basic grid data structures.
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision phys_trac_mask
integer, parameter rpxmin
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
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.
double precision unit_temperature
Physical scaling factor for temperature.
integer, parameter rpxmax
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...
procedure(sub_get_v), pointer phys_get_v
subroutine trace_field_single(xf, wP, wL, dL, numP, nwP, nwL, forward, ftype, tcondi)
subroutine get_k(xfn, igrid, K, ixIL, dxbD, ftype)
subroutine find_next_grid(igrid, igrid_next, ipe_next, xf1, newpe, stopT)
subroutine get_t_loc_trac(igrid, xloc, Tloc, ixIL, dxbD)
subroutine trace_field_multi(xfm, wPm, wLm, dL, numL, numP, nwP, nwL, forwardm, ftype, tcondi)
subroutine find_points_in_pe(igrid, ipoint_in, xf, wP, wL, dL, numP, nwP, nwL, forward, ftype, tcondi, statusF)
subroutine find_points_interp(igrid, ip_in, ip_out, xf, wP, wL, numP, nwP, nwL, dL, forward, ftype, tcondi)
Module with all the methods that users can customize in AMRVAC.
procedure(set_field_w), pointer usr_set_field_w
procedure(set_field), pointer usr_set_field