9 double precision :: dL,Tmax,trac_delta,T_bott
10 double precision,
allocatable :: xFi(:,:)
13 double precision :: dxT^D
14 double precision :: xTmin(ndim),xTmax(ndim)
15 double precision,
allocatable :: xT(:^D&,:)
17 integer,
allocatable :: trac_grid(:),ground_grid(:)
18 integer :: ngrid_trac,ngrid_ground
19 logical,
allocatable :: trac_pe(:)
23 subroutine init_trac_line(mask)
24 logical,
intent(in) :: mask
25 integer :: refine_factor,ix^D,ix(ndim),j,iFL,numL(ndim),finegrid
26 double precision :: lengthFL
27 double precision :: xprobmin(ndim),xprobmax(ndim),domain_nx(ndim)
33 ^d&xprobmin(^d)=xprobmin^d\
34 ^d&xprobmax(^d)=xprobmax^d\
35 ^d&domain_nx(^d)=domain_nx^d\
36 dl=(xprobmax(ndim)-xprobmin(ndim))/(domain_nx(ndim)*refine_factor)
46 numlp=floor(lengthfl/dl)
52 numl(j)=floor((xprobmax(j)-xprobmin(j))/dl/finegrid)
55 allocate(xfi(numfl,ndim))
56 xfi(:,ndim)=xprobmin(ndim)+dl/50.d0
57 {
do ix^db=1,numl(^db)\}
61 ifl=ifl+(ix(j)-(ndim-1-j))*(numl(j))**(ndim-1-j)
63 xfi(ifl,1:ndim-1)=xprobmin(1:ndim-1)+finegrid*dl*ix(1:ndim-1)-finegrid*dl/2.d0
66 if(
mype .eq. 0)
write(*,*)
'NOTE: 2D TRAC method take the y-dir == grav-dir'
69 if(
mype .eq. 0)
write(*,*)
'NOTE: 3D TRAC method take the z-dir == grav-dir'
72 memxfi=floor(8*numfl*numlp*ndim/1e6)
73 if (
mype==0)
write(*,*)
'Memory requirement for each processor in TRAC:'
74 if (
mype==0)
write(*,*) memxfi,
' MB'
78 end subroutine init_trac_line
80 subroutine init_trac_block(mask)
81 logical,
intent(in) :: mask
82 integer :: refine_factor,finegrid,iFL,j
84 integer :: numL(ndim),ix(ndim)
85 double precision :: lengthFL
86 double precision :: ration,a0
87 double precision :: xprobmin(ndim),xprobmax(ndim),dxT(ndim)
89 refine_factor=2**(refine_max_level-1)
90 ^d&xprobmin(^d)=xprobmin^d\
91 ^d&xprobmax(^d)=xprobmax^d\
92 ^d&dxt^d=(xprobmax^d-xprobmin^d)/(domain_nx^d*refine_factor/block_nx^d)\
94 finegrid=phys_trac_finegrid
99 dl=min(dxt^d)/finegrid
102 ^d&xtmin(^d)=xprobmin^d\
103 ^d&xtmax(^d)=xprobmax^d\
104 if(mask) xtmax(ndim)=phys_trac_mask
107 lengthfl=maxval(xtmax-xprobmin)*3.d0
109 lengthfl=maxval(xprobmax-xprobmin)*3.d0
111 numlp=floor(lengthfl/dl)
112 ^d&numxt^d=ceiling((xtmax(^d)-xtmin(^d)-smalldouble)/dxt^d)\
113 allocate(xt(numxt^d,ndim))
117 xt(j^d%ixT^s,^d)=(j-0.5d0)*dxt^d+xtmin(^d)
119 if(mask) xtmax(ndim)=maxval(xt(:^d&,ndim))+half*dxt(ndim)
124 numl(j)=floor((xprobmax(j)-xprobmin(j))/dl)
127 allocate(xfi(numfl,ndim))
128 xfi(:,ndim)=xprobmin(ndim)+dl/50.d0
129 {
do ix^db=1,numl(^db)\}
133 ifl=ifl+(ix(j)-(ndim-1-j))*(numl(j))**(ndim-1-j)
135 xfi(ifl,1:ndim-1)=xprobmin(1:ndim-1)+dl*ix(1:ndim-1)-dl/2.d0
138 if(mype .eq. 0)
write(*,*)
'NOTE: 2D TRAC method take the y-dir == grav-dir'
141 if(mype .eq. 0)
write(*,*)
'NOTE: 3D TRAC method take the z-dir == grav-dir'
143 end subroutine init_trac_block
145 subroutine trac_simple(tco_global,trac_alfa,T_peak)
146 double precision,
intent(in) :: tco_global, trac_alfa,T_peak
147 integer :: iigrid, igrid
149 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
151 ps(igrid)%special_values(1)=tco_global
153 if(ps(igrid)%special_values(1)<trac_alfa*ps(igrid)%special_values(2))
then
154 ps(igrid)%special_values(1)=trac_alfa*ps(igrid)%special_values(2)
156 if(ps(igrid)%special_values(1) .lt. t_bott)
then
157 ps(igrid)%special_values(1)=t_bott
158 else if(ps(igrid)%special_values(1) .gt. 0.2d0*t_peak)
then
159 ps(igrid)%special_values(1)=0.2d0*t_peak
161 ps(igrid)%wextra(ixg^t,iw_tcoff)=ps(igrid)%special_values(1)
163 ps(igrid)%special_values(2)=ps(igrid)%special_values(1)
165 end subroutine trac_simple
167 subroutine ltrac(T_peak)
168 double precision,
intent(in) :: T_peak
169 integer :: iigrid, igrid
170 integer :: ixO^L,trac_tcoff
174 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
175 where(ps(igrid)%wextra(ixo^s,trac_tcoff) .lt. t_bott)
176 ps(igrid)%wextra(ixo^s,trac_tcoff)=t_bott
177 else where(ps(igrid)%wextra(ixo^s,trac_tcoff) .gt. 0.2d0*t_peak)
178 ps(igrid)%wextra(ixo^s,trac_tcoff)=0.2d0*t_peak
183 subroutine tracl(mask,T_peak)
184 logical,
intent(in) :: mask
185 double precision,
intent(in) :: T_peak
187 integer :: iigrid, igrid
188 double precision :: xF(numFL,numLP,ndim)
189 integer :: numR(numFL),ix^L
190 double precision :: Tlcoff(numFL)
191 integer :: ipel(numFL,numLP),igridl(numFL,numLP)
192 logical :: forwardl(numFL)
201 call mpi_barrier(icomm,ierrmpi)
208 call get_btracing_dir(ipel,igridl,forwardl)
210 call get_tcoff_line(xf,numr,tlcoff,ipel,igridl,forwardl,mask)
212 call init_trac_tcoff()
214 call interp_tcoff(xf,ipel,igridl,numr,tlcoff)
216 call mpi_barrier(icomm,ierrmpi)
221 subroutine tracb(mask,T_peak)
222 logical,
intent(in) :: mask
223 double precision,
intent(in) :: T_peak
224 integer :: peArr(numxT^D),gdArr(numxT^D),numR(numFL)
225 double precision :: Tcoff(numxT^D),Tcmax(numxT^D),Bdir(numxT^D,ndim)
226 double precision :: xF(numFL,numLP,ndim),Tcoff_line(numFL)
227 integer :: xpe(numFL,numLP,2**ndim)
228 integer :: xgd(numFL,numLP,2**ndim)
236 call block_estable(mask,tcoff,tcmax,bdir,pearr,gdarr)
240 call block_trace_mfl(mask,tcoff,tcoff_line,tcmax,bdir,pearr,gdarr,xf,numr,xpe,xgd)
241 call block_interp_grid(mask,xf,numr,xpe,xgd,tcoff_line)
244 subroutine block_estable(mask,Tcoff,Tcmax,Bdir,peArr,gdArr)
246 double precision :: Tcoff(numxT^D),Tcoff_recv(numxT^D)
247 double precision :: Tcmax(numxT^D),Tcmax_recv(numxT^D)
248 double precision :: Bdir(numxT^D,ndim),Bdir_recv(numxT^D,ndim)
249 integer :: peArr(numxT^D),peArr_recv(numxT^D)
250 integer :: gdArr(numxT^D),gdArr_recv(numxT^D)
251 integer :: xc^L,xd^L,ix^D
252 integer :: iigrid,igrid,numxT,intab
253 double precision :: xb^L
261 xcmin^d=nghostcells+1\
262 xcmax^d=block_nx^d+nghostcells\
263 do iigrid=1,igridstail; igrid=igrids(iigrid);
264 ps(igrid)%wextra(:^d&,iw_tweight)=zero
265 ps(igrid)%wextra(:^d&,iw_tcoff)=zero
266 ^d&xbmin^d=rnode(rpxmin^d_,igrid)-xtmin(^d)\
267 ^d&xbmax^d=rnode(rpxmax^d_,igrid)-xtmin(^d)\
268 xdmin^d=nint(xbmin^d/dxt^d)+1\
269 xdmax^d=ceiling((xbmax^d-smalldouble)/dxt^d)\
270 {
do ix^d=xdmin^d,xdmax^d \}
272 {
if (ix^d .le. numxt^d) intab=intab+1 \}
273 if(intab .eq. ndim)
then
275 tcoff(ix^d)=max(tcoff(ix^d),ps(igrid)%special_values(1))
276 tcmax(ix^d)=ps(igrid)%special_values(2)
278 bdir(ix^d,1:ndim)=ps(igrid)%special_values(3:3+ndim-1)+2.d0
284 call mpi_barrier(icomm,ierrmpi)
286 call mpi_allreduce(pearr,pearr_recv,numxt,mpi_integer,&
287 mpi_max,icomm,ierrmpi)
288 call mpi_allreduce(gdarr,gdarr_recv,numxt,mpi_integer,&
289 mpi_max,icomm,ierrmpi)
290 call mpi_allreduce(tcoff,tcoff_recv,numxt,mpi_double_precision,&
291 mpi_max,icomm,ierrmpi)
292 call mpi_allreduce(bdir,bdir_recv,numxt*ndim,mpi_double_precision,&
293 mpi_max,icomm,ierrmpi)
295 call mpi_allreduce(tcmax,tcmax_recv,numxt,mpi_double_precision,&
296 mpi_max,icomm,ierrmpi)
302 if(.not. mask) tcmax=tcmax_recv
303 end subroutine block_estable
305 subroutine block_trace_mfl(mask,Tcoff,Tcoff_line,Tcmax,Bdir,peArr,gdArr,xF,numR,xpe,xgd)
306 integer :: i,j,k,k^D,ix_next^D
307 logical :: mask,flag,first
308 double precision :: Tcoff(numxT^D),Tcoff_line(numFL)
309 double precision :: Tcmax(numxT^D),Tcmax_line(numFL)
310 double precision :: xF(numFL,numLP,ndim)
311 integer :: ix_mod(ndim,2),numR(numFL)
312 double precision :: alfa_mod(ndim,2)
313 double precision :: nowpoint(ndim),nowgridc(ndim)
314 double precision :: Bdir(numxT^D,ndim)
315 double precision :: init_dir,now_dir1(ndim),now_dir2(ndim)
316 integer :: peArr(numxT^D),xpe(numFL,numLP,2**ndim)
317 integer :: gdArr(numxT^D),xgd(numFL,numLP,2**ndim)
321 ^d&k^d=ceiling((xfi(i,^d)-xtmin(^d)-smalldouble)/dxt^d)\
322 tcoff_line(i)=tcoff(k^d)
323 if(.not. mask) tcmax_line(i)=tcmax(k^d)
328 nowpoint(:)=xf(i,j,:)
329 nowgridc(:)=xt(ix_next^d,:)
332 call rk_bdir(nowgridc,nowpoint,ix_next^d,now_dir1,bdir,&
333 ix_mod,first,init_dir)
335 call rk_bdir(nowgridc,nowpoint,ix_next^d,now_dir1,bdir,&
339 xgd(i,j,1)=gdarr(ix_mod(1,1),ix_mod(2,1))
340 xgd(i,j,2)=gdarr(ix_mod(1,2),ix_mod(2,1))
341 xgd(i,j,3)=gdarr(ix_mod(1,1),ix_mod(2,2))
342 xgd(i,j,4)=gdarr(ix_mod(1,2),ix_mod(2,2))
343 xpe(i,j,1)=pearr(ix_mod(1,1),ix_mod(2,1))
344 xpe(i,j,2)=pearr(ix_mod(1,2),ix_mod(2,1))
345 xpe(i,j,3)=pearr(ix_mod(1,1),ix_mod(2,2))
346 xpe(i,j,4)=pearr(ix_mod(1,2),ix_mod(2,2))
349 xgd(i,j,1)=gdarr(ix_mod(1,1),ix_mod(2,1),ix_mod(3,1))
350 xgd(i,j,2)=gdarr(ix_mod(1,2),ix_mod(2,1),ix_mod(3,1))
351 xgd(i,j,3)=gdarr(ix_mod(1,1),ix_mod(2,2),ix_mod(3,1))
352 xgd(i,j,4)=gdarr(ix_mod(1,2),ix_mod(2,2),ix_mod(3,1))
353 xgd(i,j,5)=gdarr(ix_mod(1,1),ix_mod(2,1),ix_mod(3,2))
354 xgd(i,j,6)=gdarr(ix_mod(1,2),ix_mod(2,1),ix_mod(3,2))
355 xgd(i,j,7)=gdarr(ix_mod(1,1),ix_mod(2,2),ix_mod(3,2))
356 xgd(i,j,8)=gdarr(ix_mod(1,2),ix_mod(2,2),ix_mod(3,2))
357 xpe(i,j,1)=pearr(ix_mod(1,1),ix_mod(2,1),ix_mod(3,1))
358 xpe(i,j,2)=pearr(ix_mod(1,2),ix_mod(2,1),ix_mod(3,1))
359 xpe(i,j,3)=pearr(ix_mod(1,1),ix_mod(2,2),ix_mod(3,1))
360 xpe(i,j,4)=pearr(ix_mod(1,2),ix_mod(2,2),ix_mod(3,1))
361 xpe(i,j,5)=pearr(ix_mod(1,1),ix_mod(2,1),ix_mod(3,2))
362 xpe(i,j,6)=pearr(ix_mod(1,2),ix_mod(2,1),ix_mod(3,2))
363 xpe(i,j,7)=pearr(ix_mod(1,1),ix_mod(2,2),ix_mod(3,2))
364 xpe(i,j,8)=pearr(ix_mod(1,2),ix_mod(2,2),ix_mod(3,2))
366 nowpoint(:)=nowpoint(:)+init_dir*now_dir1*dl
367 {
if(nowpoint(^d) .gt. xtmax(^d) .or. nowpoint(^d) .lt. xtmin(^d))
then
370 if(mask .and. nowpoint(ndim) .gt. phys_trac_mask)
then
375 ^d&ix_next^d=ceiling((nowpoint(^d)-xtmin(^d)-smalldouble)/dxt^d)\
376 nowgridc(:)=xt(ix_next^d,:)
377 call rk_bdir(nowgridc,nowpoint,ix_next^d,now_dir2,bdir,&
379 xf(i,j+1,:)=xf(i,j,:)+init_dir*dl*half*(now_dir1+now_dir2)
380 {
if(xf(i,j+1,^d) .gt. xtmax(^d) .or. xf(i,j+1,^d) .lt. xtmin(^d))
then
383 if(mask .and. xf(i,j+1,ndim) .gt. phys_trac_mask)
then
387 ^d&ix_next^d=ceiling((xf(i,j+1,^d)-xtmin(^d)-smalldouble)/dxt^d)\
389 tcoff_line(i)=max(tcoff_line(i),tcoff(ix_next^d))
390 if(.not.mask) tcmax_line(i)=max(tcmax_line(i),tcmax(ix_next^d))
396 if(tcoff_line(i) .gt. tmax*0.2d0)
then
397 tcoff_line(i)=tmax*0.2d0
400 if(tcoff_line(i) .gt. tcmax_line(i)*0.2d0)
then
401 tcoff_line(i)=tcmax_line(i)*0.2d0
405 end subroutine block_trace_mfl
407 subroutine rk_bdir(nowgridc,nowpoint,ix_next^D,now_dir,Bdir,ix_mod,first,init_dir)
408 double precision :: nowpoint(ndim),nowgridc(ndim)
409 integer :: ix_mod(ndim,2)
410 double precision :: alfa_mod(ndim,2)
411 integer :: ix_next^D,k^D
412 double precision :: now_dir(ndim)
413 double precision :: Bdir(numxT^D,ndim)
415 double precision,
optional :: init_dir
417 {
if(nowpoint(^d) .gt. xtmin(^d)+half*dxt^d .and. nowpoint(^d) .lt. xtmax(^d)-half*dxt^d)
then
418 if(nowpoint(^d) .le. nowgridc(^d))
then
419 ix_mod(^d,1)=ix_next^d-1
420 ix_mod(^d,2)=ix_next^d
421 alfa_mod(^d,1)=abs(nowgridc(^d)-nowpoint(^d))/dxt^d
422 alfa_mod(^d,2)=one-alfa_mod(^d,1)
424 ix_mod(^d,1)=ix_next^d
425 ix_mod(^d,2)=ix_next^d+1
426 alfa_mod(^d,2)=abs(nowgridc(^d)-nowpoint(^d))/dxt^d
427 alfa_mod(^d,1)=one-alfa_mod(^d,2)
430 ix_mod(^d,:)=ix_next^d
437 now_dir=now_dir + bdir(ix_mod(1,k1),ix_mod(2,k2),:)*alfa_mod(1,k1)*alfa_mod(2,k2)
445 now_dir=now_dir + bdir(ix_mod(1,k1),ix_mod(2,k2),ix_mod(3,k3),:)&
446 *alfa_mod(1,k1)*alfa_mod(2,k2)*alfa_mod(3,k3)
451 if(
present(init_dir))
then
452 init_dir=sign(one,now_dir(ndim))
454 end subroutine rk_bdir
456 subroutine block_interp_grid(mask,xF,numR,xpe,xgd,Tcoff_line)
458 double precision :: xF(numFL,numLP,ndim)
459 integer :: numR(numFL)
460 integer :: xpe(numFL,numLP,2**ndim)
461 integer :: xgd(numFL,numLP,2**ndim)
462 double precision :: Tcoff_line(numFL)
463 double precision :: weightIndex,weight,ds
464 integer :: i,j,k,igrid,iigrid,ixO^L,ixc^L,ixc^D
465 double precision :: dxMax^D,dxb^D
474 if(mype .eq. xpe(i,j,k))
then
476 if(igrid .le. igrids(igridstail))
then
477 ^d&dxb^d=rnode(rpdx^d_,igrid)\
478 ^d&ixcmin^d=floor((xf(i,j,^d)-dxmax^d-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d\
479 ^d&ixcmax^d=floor((xf(i,j,^d)+dxmax^d-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d\
480 {
if (ixcmin^d<ixomin^d) ixcmin^d=ixomin^d\}
481 {
if (ixcmax^d>ixomax^d) ixcmax^d=ixomax^d\}
482 {
do ixc^d=ixcmin^d,ixcmax^d\}
484 {ds=ds+(xf(i,j,^d)-ps(igrid)%x(ixc^dd,^d))**2\}
486 if(ds .le. 0.099d0*dl)
then
487 weight=(1/(0.099d0*dl))**weightindex
489 weight=(1/ds)**weightindex
491 ps(igrid)%wextra(ixc^d,iw_tweight)=ps(igrid)%wextra(ixc^d,iw_tweight)+weight
492 ps(igrid)%wextra(ixc^d,iw_tcoff)=ps(igrid)%wextra(ixc^d,iw_tcoff)+weight*tcoff_line(i)
495 call mpistop(
"we need to check here 366Line in mod_trac.t")
502 do iigrid=1,igridstail; igrid=igrids(iigrid);
503 where (ps(igrid)%wextra(ixo^s,iw_tweight)>0.d0)
504 ps(igrid)%wextra(ixo^s,iw_tcoff)=ps(igrid)%wextra(ixo^s,iw_tcoff)/ps(igrid)%wextra(ixo^s,iw_tweight)
506 ps(igrid)%wextra(ixo^s,iw_tcoff)=0.2d0*tmax
509 end subroutine block_interp_grid
524 subroutine init_trac_tcoff()
525 integer :: ixI^L,ixO^L,igrid,iigrid
530 do iigrid=1,igridstail; igrid=igrids(iigrid);
531 ps(igrid)%wextra(ixi^s,iw_tcoff)=0.d0
532 ps(igrid)%wextra(ixi^s,iw_tweight)=0.d0
534 end subroutine init_trac_tcoff
536 subroutine update_pegrid()
545 call traverse_gridtable()
547 end subroutine update_pegrid
549 subroutine traverse_gridtable()
552 double precision :: dxb^D,xb^L
553 integer :: iigrid,igrid,j
554 logical,
allocatable :: trac_pe_recv(:)
555 double precision :: hcmax_bt
557 allocate(trac_pe_recv(
npe))
561 do iigrid=1,igridstail; igrid=igrids(iigrid);
565 ngrid_trac=ngrid_trac+1
566 trac_grid(ngrid_trac)=igrid
567 if (xbmin^nd<hcmax_bt)
then
568 ngrid_ground=ngrid_ground+1
569 ground_grid(ngrid_ground)=igrid
574 call mpi_allreduce(trac_pe,trac_pe_recv,
npe,mpi_logical,mpi_lor,
icomm,
ierrmpi)
577 deallocate(trac_pe_recv)
578 end subroutine traverse_gridtable
580 subroutine get_te_grid()
581 integer :: ixI^L,ixO^L,igrid,iigrid,j
582 double precision :: rho(ixM^T)
588 call phys_get_pthermal(ps(igrid)%w,ps(igrid)%x,ixi^l,ixi^l,ps(igrid)%wextra(ixi^s,iw_tcoff))
589 call phys_get_rho(ps(igrid)%w,ps(igrid)%x,ixi^l,ixi^l,rho)
590 ps(igrid)%wextra(ixi^s,iw_tcoff)=ps(igrid)%wextra(ixi^s,iw_tcoff)/rho(ixi^s)
592 end subroutine get_te_grid
594 subroutine get_btracing_dir(ipel,igridl,forwardl)
595 integer :: ipel(numFL,numLP),igridl(numFL,numLP)
596 logical :: forwardl(numFL)
597 integer :: igrid,ixO^L,iFL,j,ix^D,idir,ixb^D,ixbb^D
598 double precision :: xb^L,dxb^D,xd^D,factor,Bh
599 integer :: numL(ndim),ixmin(ndim),ixmax(ndim),ix(ndim)
600 logical :: forwardRC(numFL)
610 ^d&dxb^d=rnode(rpdx^d_,igrid);
611 ^d&xbmin^d=rnode(rpxmin^d_,igrid);
612 ^d&xbmax^d=rnode(rpxmax^d_,igrid);
613 ^d&ixmin(^d)=floor((xbmin^d-xprobmin^d)/(phys_trac_finegrid*dl))+1;
614 ^d&ixmax(^d)=floor((xbmax^d-xprobmin^d)/(phys_trac_finegrid*dl));
615 ^d&numl(^d)=floor((xprobmax^d-xprobmin^d)/(phys_trac_finegrid*dl));
619 {
do ix^db=ixmin(^db),ixmax(^db)\}
623 ifl=ifl+(ix(idir)-(ndim-1-idir))*(numl(idir))**(ndim-1-idir)
627 ^d&ixb^d=floor((xfi(ifl,^d)-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d;
628 ^d&xd^d=(xfi(ifl,^d)-ps(igrid)%x(ixb^dd,^d))/dxb^d;
631 factor={abs(1-ix^d-xd^d)*}
633 if(
allocated(iw_mag))
then
634 bh=bh+factor*(ps(igrid)%w(ixb^d+ixbb^d,iw_mag(^nd))+ps(igrid)%B0(ixb^d+ixbb^d,^nd,0))
636 bh=bh+factor*(ps(igrid)%B0(ixb^d+ixbb^d,^nd,0))
639 bh=bh+factor*ps(igrid)%w(ixb^d+ixbb^d,iw_mag(^nd))
645 forwardl(ifl)=.false.
649 call mpi_allreduce(forwardl,forwardrc,numfl,mpi_logical,&
650 mpi_land,icomm,ierrmpi)
652 end subroutine get_btracing_dir
654 subroutine get_tcoff_line(xFL,numR,TcoffFL,ipeFL,igridFL,forwardFL,mask)
656 double precision :: xFL(numFL,numLP,ndim)
657 integer :: numR(numFL)
658 double precision :: TcoffFL(numFL),TmaxFL(numFL)
659 integer :: ipeFL(numFL,numLP),igridFL(numFL,numLP)
660 logical :: forwardFL(numFL)
661 logical,
intent(in) :: mask
662 integer :: nwP,nwL,iFL,iLP
663 double precision :: wPm(numFL,numLP,2),wLm(numFL,1+2)
664 character(len=std_len) :: ftype,tcondi
670 call trace_field_multi(xfl,wpm,wlm,dl,numfl,numlp,nwp,nwl,forwardfl,ftype,tcondi)
672 numr(ifl)=int(wlm(ifl,1))
673 tcofffl(ifl)=wlm(ifl,2)
674 tmaxfl(ifl)=wlm(ifl,3)
676 if(tcofffl(ifl)>0.2d0*tmax) tcofffl(ifl)=0.2d0*tmax
678 tmaxfl(ifl)=wlm(ifl,3)
679 if(tcofffl(ifl)>0.2d0*tmaxfl(ifl)) tcofffl(ifl)=0.2d0*tmaxfl(ifl)
682 if(tcofffl(ifl)<t_bott) tcofffl(ifl)=t_bott
686 if (numr(ifl)>0)
then
688 ipefl(ifl,ilp)=int(wpm(ifl,ilp,1))
689 igridfl(ifl,ilp)=int(wpm(ifl,ilp,2))
693 end subroutine get_tcoff_line
695 subroutine interp_tcoff(xF,ipel,igridl,numR,Tlcoff)
696 double precision :: xF(numFL,numLP,ndim)
697 integer :: numR(numFL),ipel(numFL,numLP),igridl(numFL,numLP)
698 double precision :: Tlcoff(numFL)
699 integer :: iFL,iLP,ixO^L,ixI^L,ixc^L,ixb^L,ixc^D
700 integer :: igrid,j,ipmin,ipmax,igrid_nb
701 double precision :: dxb^D,dxMax^D,xb^L,Tcnow
702 double precision :: xFnow(ndim),xc(ndim)
703 integer :: weightIndex,idn^D,ixmax^ND
704 double precision :: ds,weight
716 do while (ilp<=numr(ifl))
719 do while (ipel(ifl,ipmin)/=mype .and. ipmin<=numr(ifl))
722 igrid=igridl(ifl,ipmin)
724 do while (ipel(ifl,ipmax)==mype .and. igridl(ifl,ipmax+1)==igrid .and. ipmax<numr(ifl))
729 ^d&dxb^d=rnode(rpdx^d_,igrid);
733 xfnow(:)=xf(ifl,ilp,:)
734 ^d&ixbmin^d=floor((xfnow(^d)-dxmax^d-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d;
735 ^d&ixbmax^d=floor((xfnow(^d)+dxmax^d-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d;
738 {ixcmin^d=max(ixbmin^d,ixomin^d)\}
739 {ixcmax^d=min(ixbmax^d,ixomax^d)\}
740 xbmin^nd=rnode(rpxmin^nd_,igrid)
741 xbmax^nd=rnode(rpxmax^nd_,igrid)
742 ixmax^nd=floor((phys_trac_mask-xbmin^nd)/dxb^nd)+ixomin^nd
743 if (xbmax^nd>phys_trac_mask) ixcmax^nd=min(ixmax^nd,ixcmax^nd)
744 {
do ixc^d=ixcmin^d,ixcmax^d\}
746 {ds=ds+(xfnow(^d)-ps(igrid)%x(ixc^dd,^d))**2\}
748 if(ds<1.0d-2*dxb1)
then
749 weight=(1/(1.0d-2*dxb1))**weightindex
751 weight=(1/ds)**weightindex
753 ps(igrid)%wextra(ixc^d,iw_tweight)=ps(igrid)%wextra(ixc^d,iw_tweight)+weight
754 ps(igrid)%wextra(ixc^d,iw_tcoff)=ps(igrid)%wextra(ixc^d,iw_tcoff)+weight*tcnow
759 if (ixbmin^d<ixomin^d)
then
762 if (neighbor(2,idn^dd,igrid)==mype .and. neighbor_type(idn^dd,igrid)==neighbor_sibling)
then
763 igrid_nb=neighbor(1,idn^dd,igrid)
764 ixcmin^dd=max(ixbmin^dd,ixomin^dd);
765 ixcmax^dd=min(ixbmax^dd,ixomax^dd);
766 ixcmin^d=ixomax^d+(ixbmin^d-ixomin^d)
768 xbmin^nd=rnode(rpxmin^nd_,igrid_nb)
769 xbmax^nd=rnode(rpxmax^nd_,igrid_nb)
770 ixmax^nd=floor((phys_trac_mask-xbmin^nd)/dxb^nd)+ixomin^nd
771 if (xbmax^nd>phys_trac_mask) ixcmax^nd=min(ixmax^nd,ixcmax^nd)
773 {
do ixc^dd=ixcmin^dd,ixcmax^dd;}
775 xc(:)=ps(igrid_nb)%x(ixc^dd,:)
776 {ds=ds+(xfnow(^dd)-xc(^dd))**2;}
778 if(ds<1.0d-2*dxb1)
then
779 weight=(1/(1.0d-2*dxb1))**weightindex
781 weight=(1/ds)**weightindex
783 ps(igrid_nb)%wextra(ixc^dd,iw_tweight)=ps(igrid_nb)%wextra(ixc^dd,iw_tweight)+weight
784 ps(igrid_nb)%wextra(ixc^dd,iw_tcoff)=ps(igrid_nb)%wextra(ixc^dd,iw_tcoff)+weight*tcnow
789 if (ixbmax^d>ixomin^d)
then
792 if (neighbor(2,idn^dd,igrid)==mype .and. neighbor_type(idn^dd,igrid)==neighbor_sibling)
then
793 igrid_nb=neighbor(1,idn^dd,igrid)
794 xbmin^nd=rnode(rpxmin^nd_,igrid_nb)
795 if (xbmin^nd<phys_trac_mask)
then
796 ixcmin^dd=max(ixbmin^dd,ixomin^dd);
797 ixcmax^dd=min(ixbmax^dd,ixomax^dd);
799 ixcmax^d=ixomin^d+(ixbmax^d-ixomax^d)
800 xbmax^nd=rnode(rpxmax^nd_,igrid_nb)
801 ixmax^nd=floor((phys_trac_mask-xbmin^nd)/dxb^nd)+ixomin^nd
802 if (xbmax^nd>phys_trac_mask) ixcmax^nd=min(ixmax^nd,ixcmax^nd)
804 {
do ixc^dd=ixcmin^dd,ixcmax^dd;}
806 xc(:)=ps(igrid_nb)%x(ixc^dd,:)
807 {ds=ds+(xfnow(^dd)-xc(^dd))**2;}
809 if(ds<1.0d-2*dxb1)
then
810 weight=(1/(1.0d-2*dxb1))**weightindex
812 weight=(1/ds)**weightindex
814 ps(igrid_nb)%wextra(ixc^dd,iw_tweight)=ps(igrid_nb)%wextra(ixc^dd,iw_tweight)+weight
815 ps(igrid_nb)%wextra(ixc^dd,iw_tcoff)=ps(igrid_nb)%wextra(ixc^dd,iw_tcoff)+weight*tcnow
829 where(ps(igrid)%wextra(ixo^s,iw_tweight)>0.d0)
830 ps(igrid)%wextra(ixo^s,iw_tcoff)=ps(igrid)%wextra(ixo^s,iw_tcoff)/ps(igrid)%wextra(ixo^s,iw_tweight)
832 ps(igrid)%wextra(ixo^s,iw_tcoff)=t_bott
835 end subroutine interp_tcoff
837 subroutine trac_after_setdt(tco,trac_alfa,T_peak, T_bott_in)
838 double precision,
intent(in) :: trac_alfa,tco,T_peak, T_bott_in
841 select case(phys_trac_type)
847 call trac_simple(tco,trac_alfa,t_peak)
854 call tracl(.false.,t_peak)
857 call tracb(.false.,t_peak)
860 call tracl(.true.,t_peak)
863 call tracb(.true.,t_peak)
865 call mpistop(
"undefined TRAC method type")
867 end subroutine trac_after_setdt
876 if(
mype .eq. 0)
write(*,*)
'Using TRACL(ine) global method'
877 if(
mype .eq. 0)
write(*,*)
'By default, magnetic field lines are traced every 4 grid cells'
878 call init_trac_line(.false.)
881 if(
mype .eq. 0)
write(*,*)
'Using TRACB(lock) global method'
882 if(
mype .eq. 0)
write(*,*)
'Currently, only valid in Cartesian uniform settings'
883 if(
mype .eq. 0)
write(*,*)
'By default, magnetic field lines are traced every 4 grid cells'
884 call init_trac_block(.false.)
887 if(
mype .eq. 0)
write(*,*)
'Using TRACL(ine) method with a mask'
888 if(
mype .eq. 0)
write(*,*)
'By default, magnetic field lines are traced every 4 grid cells'
889 call init_trac_line(.true.)
892 if(
mype .eq. 0)
write(*,*)
'Using TRACB(lock) method with a mask'
893 if(
mype .eq. 0)
write(*,*)
'Currently, only valid in Cartesian uniform settings'
894 if(
mype .eq. 0)
write(*,*)
'By default, magnetic field lines are traced every 4 grid cells'
895 call init_trac_block(.true.)
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 domain_nx
number of cells for each dimension in level-one mesh
double precision phys_trac_mask
integer, parameter rpxmin
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
logical phys_trac
Use TRAC for MHD or 1D HD.
integer refine_max_level
Maximal number of AMR levels.
integer max_blocks
The maximum number of grid blocks in a processor.
integer phys_trac_finegrid
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_trac_after_setdt), pointer phys_trac_after_setdt
subroutine, public initialize_trac_after_settree
subroutine trace_field_multi(xfm, wpm, wlm, dl, numl, nump, nwp, nwl, forwardm, ftype, tcondi)