MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_trac.t
Go to the documentation of this file.
1 module mod_trac
3  use mod_mhd
4  implicit none
5  ! common
6  integer :: numfl,numlp
7  double precision :: dl,tmax,trac_delta,t_bott
8  double precision, allocatable :: xfi(:,:)
9  ! TRACB
10  integer :: numxt^d
11  double precision :: dxt^d
12  double precision :: xtmin(ndim),xtmax(ndim)
13  double precision, allocatable :: xt(:^d&,:)
14  ! TRACL mask
15  integer, allocatable :: trac_grid(:),ground_grid(:)
17  logical, allocatable :: trac_pe(:)
18 contains
19 
20  subroutine init_trac_line(mask)
21  logical, intent(in) :: mask
22  integer :: refine_factor,ix^D,ix(ndim),j,iFL,numL(ndim),finegrid
23  double precision :: lengthFL
24  double precision :: xprobmin(ndim),xprobmax(ndim),domain_nx(ndim)
25  integer :: memxFi
26 
27  trac_delta=0.25d0
28  !----------------- init magnetic field -------------------!
29  refine_factor=2**(refine_max_level-1)
30  ^d&xprobmin(^d)=xprobmin^d\
31  ^d&xprobmax(^d)=xprobmax^d\
32  ^d&domain_nx(^d)=domain_nx^d\
33  dl=(xprobmax(ndim)-xprobmin(ndim))/(domain_nx(ndim)*refine_factor)
34  ! max length of a field line
35  if (.not. mask) phys_trac_mask=xprobmax^nd
36  {^iftwod
37  lengthfl=(phys_trac_mask-xprobmin(2))*3.d0
38  }
39  {^ifthreed
40  lengthfl=(phys_trac_mask-xprobmin(3))*3.d0
41  }
42 
43  numlp=floor(lengthfl/dl)
44  numl=1
45  numfl=1
46  do j=1,ndim-1
47  ! number of field lines, every 4 finest grid, every direction
48  finegrid=mhd_trac_finegrid
49  numl(j)=floor((xprobmax(j)-xprobmin(j))/dl/finegrid)
50  numfl=numfl*numl(j)
51  end do
52  allocate(xfi(numfl,ndim))
53  xfi(:,ndim)=xprobmin(ndim)+dl/50.d0
54  {do ix^db=1,numl(^db)\}
55  ^d&ix(^d)=ix^d\
56  ifl=0
57  do j=ndim-1,1,-1
58  ifl=ifl+(ix(j)-(ndim-1-j))*(numl(j))**(ndim-1-j)
59  end do
60  xfi(ifl,1:ndim-1)=xprobmin(1:ndim-1)+finegrid*dl*ix(1:ndim-1)-finegrid*dl/2.d0
61  {end do\}
62  {^iftwod
63  if(mype .eq. 0) write(*,*) 'NOTE: 2D TRAC method take the y-dir == grav-dir'
64  }
65  {^ifthreed
66  if(mype .eq. 0) write(*,*) 'NOTE: 3D TRAC method take the z-dir == grav-dir'
67  }
68 
69  memxfi=floor(8*numfl*numlp*ndim/1e6)
70  if (mype==0) write(*,*) 'Memory requirement for each processor in TRAC:'
71  if (mype==0) write(*,*) memxfi,' MB'
72 
74 
75  end subroutine init_trac_line
76 
77  subroutine init_trac_block(mask)
78  logical, intent(in) :: mask
79  integer :: refine_factor,finegrid,iFL,j
80  integer :: ix^D,ixT^L
81  integer :: numL(ndim),ix(ndim)
82  double precision :: lengthFL
83  double precision :: ration,a0
84  double precision :: xprobmin(ndim),xprobmax(ndim),dxT(ndim)
85 
86  refine_factor=2**(refine_max_level-1)
87  ^d&xprobmin(^d)=xprobmin^d\
88  ^d&xprobmax(^d)=xprobmax^d\
89  ^d&dxt^d=(xprobmax^d-xprobmin^d)/(domain_nx^d*refine_factor/block_nx^d)\
90  ^d&dxt(ndim)=dxt^d\
91  finegrid=mhd_trac_finegrid
92  {^ifoned
93  dl=dxt1/finegrid
94  }
95  {^nooned
96  dl=min(dxt^d)/finegrid
97  }
98  ! table for interpolation
99  ^d&xtmin(^d)=xprobmin^d\
100  ^d&xtmax(^d)=xprobmax^d\
101  if(mask) xtmax(ndim)=phys_trac_mask
102  ! max length of a field line
103  if(mask) then
104  lengthfl=maxval(xtmax-xprobmin)*3.d0
105  else
106  lengthfl=maxval(xprobmax-xprobmin)*3.d0
107  end if
108  numlp=floor(lengthfl/dl)
109  ^d&numxt^d=ceiling((xtmax(^d)-xtmin(^d)-smalldouble)/dxt^d)\
110  allocate(xt(numxt^d,ndim))
111  ixtmin^d=1\
112  ixtmax^d=numxt^d\
113  {do j=1,numxt^d
114  xt(j^d%ixT^s,^d)=(j-0.5d0)*dxt^d+xtmin(^d)
115  end do\}
116  if(mask) xtmax(ndim)=maxval(xt(:^d&,ndim))+half*dxt(ndim)
117  numl=1
118  numfl=1
119  do j=1,ndim-1
120  ! number of field lines, every 4 finest grid, every direction
121  numl(j)=floor((xprobmax(j)-xprobmin(j))/dl)
122  numfl=numfl*numl(j)
123  end do
124  allocate(xfi(numfl,ndim))
125  xfi(:,ndim)=xprobmin(ndim)+dl/50.d0
126  {do ix^db=1,numl(^db)\}
127  ^d&ix(^d)=ix^d\
128  ifl=0
129  do j=ndim-1,1,-1
130  ifl=ifl+(ix(j)-(ndim-1-j))*(numl(j))**(ndim-1-j)
131  end do
132  xfi(ifl,1:ndim-1)=xprobmin(1:ndim-1)+dl*ix(1:ndim-1)-dl/2.d0
133  {end do\}
134  {^iftwod
135  if(mype .eq. 0) write(*,*) 'NOTE: 2D TRAC method take the y-dir == grav-dir'
136  }
137  {^ifthreed
138  if(mype .eq. 0) write(*,*) 'NOTE: 3D TRAC method take the z-dir == grav-dir'
139  }
140  end subroutine init_trac_block
141 
142  subroutine trac_simple(tco_global,trac_alfa,T_peak)
143  double precision, intent(in) :: tco_global, trac_alfa,T_peak
144  integer :: iigrid, igrid
145 
146  do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
147  {^ifoned
148  ps(igrid)%special_values(1)=tco_global
149  }
150  if(ps(igrid)%special_values(1)<trac_alfa*ps(igrid)%special_values(2)) then
151  ps(igrid)%special_values(1)=trac_alfa*ps(igrid)%special_values(2)
152  end if
153  if(ps(igrid)%special_values(1) .lt. t_bott) then
154  ps(igrid)%special_values(1)=t_bott
155  else if(ps(igrid)%special_values(1) .gt. 0.2d0*t_peak) then
156  ps(igrid)%special_values(1)=0.2d0*t_peak
157  end if
158  ps(igrid)%wextra(ixg^t,iw_tcoff)=ps(igrid)%special_values(1)
159  !> special values(2) to save old tcutoff
160  ps(igrid)%special_values(2)=ps(igrid)%special_values(1)
161  end do
162  end subroutine trac_simple
163 
164  subroutine ltrac(T_peak)
165  double precision, intent(in) :: T_peak
166  integer :: iigrid, igrid
167  integer :: ixO^L,trac_tcoff
168 
169  ixo^l=ixm^ll;
170  trac_tcoff=iw_tcoff
171  do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
172  where(ps(igrid)%wextra(ixo^s,trac_tcoff) .lt. t_bott)
173  ps(igrid)%wextra(ixo^s,trac_tcoff)=t_bott
174  else where(ps(igrid)%wextra(ixo^s,trac_tcoff) .gt. 0.2d0*t_peak)
175  ps(igrid)%wextra(ixo^s,trac_tcoff)=0.2d0*t_peak
176  end where
177  end do
178  end subroutine ltrac
179 
180  subroutine tracl(mask,T_peak)
181  logical, intent(in) :: mask
182  double precision, intent(in) :: T_peak
183  integer :: ix^D,j
184  integer :: iigrid, igrid
185  double precision :: xF(numFL,numLP,ndim)
186  integer :: numR(numFL),ix^L
187  double precision :: Tlcoff(numFL)
188  integer :: ipel(numFL,numLP),igridl(numFL,numLP)
189  logical :: forwardl(numFL)
190 
191  tmax=t_peak
192  xf=zero
193 
194  do j=1,ndim
195  xf(:,1,j)=xfi(:,j)
196  end do
197 
198  call mpi_barrier(icomm,ierrmpi)
199  ! record pe and grid for mask region
200  call update_pegrid()
201  ! get temperature for the calculation of Te gradient,
202  ! which will be stored in wextra(ixI^S,Tcoff_) temporarily
203  call get_te_grid()
204  ! find out direction for tracing B field
205  call get_btracing_dir(ipel,igridl,forwardl)
206  ! trace Bfield and get Tcoff for each field line
207  call get_tcoff_line(xf,numr,tlcoff,ipel,igridl,forwardl,mask)
208  ! init vairable Tcoff_ and Tweight_
209  call init_trac_tcoff()
210  ! get cell Tcoff via interpolation
211  call interp_tcoff(xf,ipel,igridl,numr,tlcoff)
212  ! convert primitive data back to conserved data
213  call mpi_barrier(icomm,ierrmpi)
214 
215 
216  end subroutine tracl
217 
218  subroutine tracb(mask,T_peak)
219  logical, intent(in) :: mask
220  double precision, intent(in) :: T_peak
221  integer :: peArr(numxT^D),gdArr(numxT^D),numR(numFL)
222  double precision :: Tcoff(numxT^D),Tcmax(numxT^D),Bdir(numxT^D,ndim)
223  double precision :: xF(numFL,numLP,ndim),Tcoff_line(numFL)
224  integer :: xpe(numFL,numLP,2**ndim)
225  integer :: xgd(numFL,numLP,2**ndim)
226 
227  tmax=t_peak
228  tcoff=zero
229  tcmax=zero
230  bdir=zero
231  pearr=-1
232  gdarr=-1
233  call block_estable(mask,tcoff,tcmax,bdir,pearr,gdarr)
234  xf=zero
235  numr=0
236  tcoff_line=zero
237  call block_trace_mfl(mask,tcoff,tcoff_line,tcmax,bdir,pearr,gdarr,xf,numr,xpe,xgd)
238  call block_interp_grid(mask,xf,numr,xpe,xgd,tcoff_line)
239  end subroutine tracb
240 
241  subroutine block_estable(mask,Tcoff,Tcmax,Bdir,peArr,gdArr)
242  logical :: mask
243  double precision :: Tcoff(numxT^D),Tcoff_recv(numxT^D)
244  double precision :: Tcmax(numxT^D),Tcmax_recv(numxT^D)
245  double precision :: Bdir(numxT^D,ndim),Bdir_recv(numxT^D,ndim)
246  integer :: peArr(numxT^D),peArr_recv(numxT^D)
247  integer :: gdArr(numxT^D),gdArr_recv(numxT^D)
248  integer :: xc^L,xd^L,ix^D
249  integer :: iigrid,igrid,numxT,intab
250  double precision :: xb^L
251 
252  tcoff_recv=zero
253  tcmax_recv=zero
254  bdir_recv=zero
255  pearr_recv=-1
256  gdarr_recv=-1
257  !> combine table from different processors
258  xcmin^d=nghostcells+1\
259  xcmax^d=block_nx^d+nghostcells\
260  do iigrid=1,igridstail; igrid=igrids(iigrid);
261  ps(igrid)%wextra(:^d&,tweight_)=zero
262  ps(igrid)%wextra(:^d&,tcoff_)=zero
263  ^d&xbmin^d=rnode(rpxmin^d_,igrid)-xtmin(^d)\
264  ^d&xbmax^d=rnode(rpxmax^d_,igrid)-xtmin(^d)\
265  xdmin^d=nint(xbmin^d/dxt^d)+1\
266  xdmax^d=ceiling((xbmax^d-smalldouble)/dxt^d)\
267  {do ix^d=xdmin^d,xdmax^d \}
268  intab=0
269  {if (ix^d .le. numxt^d) intab=intab+1 \}
270  if(intab .eq. ndim) then
271  !> in principle, no overlap will happen here
272  tcoff(ix^d)=max(tcoff(ix^d),ps(igrid)%special_values(1))
273  tcmax(ix^d)=ps(igrid)%special_values(2)
274  !> abs(Bdir) <= 1, so that Bdir+2 should always be positive
275  bdir(ix^d,1:ndim)=ps(igrid)%special_values(3:3+ndim-1)+2.d0
276  pearr(ix^d)=mype
277  gdarr(ix^d)=igrid
278  end if
279  {end do\}
280  end do
281  call mpi_barrier(icomm,ierrmpi)
282  numxt={numxt^d*}
283  call mpi_allreduce(pearr,pearr_recv,numxt,mpi_integer,&
284  mpi_max,icomm,ierrmpi)
285  call mpi_allreduce(gdarr,gdarr_recv,numxt,mpi_integer,&
286  mpi_max,icomm,ierrmpi)
287  call mpi_allreduce(tcoff,tcoff_recv,numxt,mpi_double_precision,&
288  mpi_max,icomm,ierrmpi)
289  call mpi_allreduce(bdir,bdir_recv,numxt*ndim,mpi_double_precision,&
290  mpi_max,icomm,ierrmpi)
291  if(.not. mask) then
292  call mpi_allreduce(tcmax,tcmax_recv,numxt,mpi_double_precision,&
293  mpi_max,icomm,ierrmpi)
294  end if
295  pearr=pearr_recv
296  gdarr=gdarr_recv
297  tcoff=tcoff_recv
298  bdir=bdir_recv-2.d0
299  if(.not. mask) tcmax=tcmax_recv
300  end subroutine block_estable
301 
302  subroutine block_trace_mfl(mask,Tcoff,Tcoff_line,Tcmax,Bdir,peArr,gdArr,xF,numR,xpe,xgd)
303  integer :: i,j,k,k^D,ix_next^D
304  logical :: mask,flag,first
305  double precision :: Tcoff(numxT^D),Tcoff_line(numFL)
306  double precision :: Tcmax(numxT^D),Tcmax_line(numFL)
307  double precision :: xF(numFL,numLP,ndim)
308  integer :: ix_mod(ndim,2),numR(numFL)
309  double precision :: alfa_mod(ndim,2)
310  double precision :: nowpoint(ndim),nowgridc(ndim)
311  double precision :: Bdir(numxT^D,ndim)
312  double precision :: init_dir,now_dir1(ndim),now_dir2(ndim)
313  integer :: peArr(numxT^D),xpe(numFL,numLP,2**ndim)
314  integer :: gdArr(numxT^D),xgd(numFL,numLP,2**ndim)
315 
316  do i=1,numfl
317  flag=.true.
318  ^d&k^d=ceiling((xfi(i,^d)-xtmin(^d)-smalldouble)/dxt^d)\
319  tcoff_line(i)=tcoff(k^d)
320  if(.not. mask) tcmax_line(i)=tcmax(k^d)
321  ix_next^d=k^d\
322  j=1
323  xf(i,j,:)=xfi(i,:)
324  do while(flag)
325  nowpoint(:)=xf(i,j,:)
326  nowgridc(:)=xt(ix_next^d,:)
327  first=.true.
328  if(j .eq. 1) then
329  call rk_bdir(nowgridc,nowpoint,ix_next^d,now_dir1,bdir,&
330  ix_mod,first,init_dir)
331  else
332  call rk_bdir(nowgridc,nowpoint,ix_next^d,now_dir1,bdir,&
333  ix_mod,first)
334  end if
335  {^iftwod
336  xgd(i,j,1)=gdarr(ix_mod(1,1),ix_mod(2,1))
337  xgd(i,j,2)=gdarr(ix_mod(1,2),ix_mod(2,1))
338  xgd(i,j,3)=gdarr(ix_mod(1,1),ix_mod(2,2))
339  xgd(i,j,4)=gdarr(ix_mod(1,2),ix_mod(2,2))
340  xpe(i,j,1)=pearr(ix_mod(1,1),ix_mod(2,1))
341  xpe(i,j,2)=pearr(ix_mod(1,2),ix_mod(2,1))
342  xpe(i,j,3)=pearr(ix_mod(1,1),ix_mod(2,2))
343  xpe(i,j,4)=pearr(ix_mod(1,2),ix_mod(2,2))
344  }
345  {^ifthreed
346  xgd(i,j,1)=gdarr(ix_mod(1,1),ix_mod(2,1),ix_mod(3,1))
347  xgd(i,j,2)=gdarr(ix_mod(1,2),ix_mod(2,1),ix_mod(3,1))
348  xgd(i,j,3)=gdarr(ix_mod(1,1),ix_mod(2,2),ix_mod(3,1))
349  xgd(i,j,4)=gdarr(ix_mod(1,2),ix_mod(2,2),ix_mod(3,1))
350  xgd(i,j,5)=gdarr(ix_mod(1,1),ix_mod(2,1),ix_mod(3,2))
351  xgd(i,j,6)=gdarr(ix_mod(1,2),ix_mod(2,1),ix_mod(3,2))
352  xgd(i,j,7)=gdarr(ix_mod(1,1),ix_mod(2,2),ix_mod(3,2))
353  xgd(i,j,8)=gdarr(ix_mod(1,2),ix_mod(2,2),ix_mod(3,2))
354  xpe(i,j,1)=pearr(ix_mod(1,1),ix_mod(2,1),ix_mod(3,1))
355  xpe(i,j,2)=pearr(ix_mod(1,2),ix_mod(2,1),ix_mod(3,1))
356  xpe(i,j,3)=pearr(ix_mod(1,1),ix_mod(2,2),ix_mod(3,1))
357  xpe(i,j,4)=pearr(ix_mod(1,2),ix_mod(2,2),ix_mod(3,1))
358  xpe(i,j,5)=pearr(ix_mod(1,1),ix_mod(2,1),ix_mod(3,2))
359  xpe(i,j,6)=pearr(ix_mod(1,2),ix_mod(2,1),ix_mod(3,2))
360  xpe(i,j,7)=pearr(ix_mod(1,1),ix_mod(2,2),ix_mod(3,2))
361  xpe(i,j,8)=pearr(ix_mod(1,2),ix_mod(2,2),ix_mod(3,2))
362  }
363  nowpoint(:)=nowpoint(:)+init_dir*now_dir1*dl
364  {if(nowpoint(^d) .gt. xtmax(^d) .or. nowpoint(^d) .lt. xtmin(^d)) then
365  flag=.false.
366  end if\}
367  if(mask .and. nowpoint(ndim) .gt. phys_trac_mask) then
368  flag=.false.
369  end if
370  if(flag) then
371  first=.false.
372  ^d&ix_next^d=ceiling((nowpoint(^d)-xtmin(^d)-smalldouble)/dxt^d)\
373  nowgridc(:)=xt(ix_next^d,:)
374  call rk_bdir(nowgridc,nowpoint,ix_next^d,now_dir2,bdir,&
375  ix_mod,first)
376  xf(i,j+1,:)=xf(i,j,:)+init_dir*dl*half*(now_dir1+now_dir2)
377  {if(xf(i,j+1,^d) .gt. xtmax(^d) .or. xf(i,j+1,^d) .lt. xtmin(^d)) then
378  flag=.false.
379  end if\}
380  if(mask .and. xf(i,j+1,ndim) .gt. phys_trac_mask) then
381  flag=.false.
382  end if
383  if(flag) then
384  ^d&ix_next^d=ceiling((xf(i,j+1,^d)-xtmin(^d)-smalldouble)/dxt^d)\
385  j=j+1
386  tcoff_line(i)=max(tcoff_line(i),tcoff(ix_next^d))
387  if(.not.mask) tcmax_line(i)=max(tcmax_line(i),tcmax(ix_next^d))
388  end if
389  end if
390  end do
391  numr(i)=j
392  if(mask) then
393  if(tcoff_line(i) .gt. tmax*0.2d0) then
394  tcoff_line(i)=tmax*0.2d0
395  end if
396  else
397  if(tcoff_line(i) .gt. tcmax_line(i)*0.2d0) then
398  tcoff_line(i)=tcmax_line(i)*0.2d0
399  end if
400  end if
401  end do
402  end subroutine block_trace_mfl
403 
404  subroutine rk_bdir(nowgridc,nowpoint,ix_next^D,now_dir,Bdir,ix_mod,first,init_dir)
405  double precision :: nowpoint(ndim),nowgridc(ndim)
406  integer :: ix_mod(ndim,2)
407  double precision :: alfa_mod(ndim,2)
408  integer :: ix_next^D,k^D
409  double precision :: now_dir(ndim)
410  double precision :: Bdir(numxT^D,ndim)
411  logical :: first
412  double precision, optional :: init_dir
413 
414  {if(nowpoint(^d) .gt. xtmin(^d)+half*dxt^d .and. nowpoint(^d) .lt. xtmax(^d)-half*dxt^d) then
415  if(nowpoint(^d) .le. nowgridc(^d)) then
416  ix_mod(^d,1)=ix_next^d-1
417  ix_mod(^d,2)=ix_next^d
418  alfa_mod(^d,1)=abs(nowgridc(^d)-nowpoint(^d))/dxt^d
419  alfa_mod(^d,2)=one-alfa_mod(^d,1)
420  else
421  ix_mod(^d,1)=ix_next^d
422  ix_mod(^d,2)=ix_next^d+1
423  alfa_mod(^d,2)=abs(nowgridc(^d)-nowpoint(^d))/dxt^d
424  alfa_mod(^d,1)=one-alfa_mod(^d,2)
425  end if
426  else
427  ix_mod(^d,:)=ix_next^d
428  alfa_mod(^d,:)=half
429  end if\}
430  now_dir=zero
431  {^iftwod
432  do k1=1,2
433  do k2=1,2
434  now_dir=now_dir + bdir(ix_mod(1,k1),ix_mod(2,k2),:)*alfa_mod(1,k1)*alfa_mod(2,k2)
435  end do
436  end do
437  }
438  {^ifthreed
439  do k1=1,2
440  do k2=1,2
441  do k3=1,2
442  now_dir=now_dir + bdir(ix_mod(1,k1),ix_mod(2,k2),ix_mod(3,k3),:)&
443  *alfa_mod(1,k1)*alfa_mod(2,k2)*alfa_mod(3,k3)
444  end do
445  end do
446  end do
447  }
448  if(present(init_dir)) then
449  init_dir=sign(one,now_dir(ndim))
450  end if
451  end subroutine rk_bdir
452 
453  subroutine block_interp_grid(mask,xF,numR,xpe,xgd,Tcoff_line)
454  logical :: mask
455  double precision :: xF(numFL,numLP,ndim)
456  integer :: numR(numFL)
457  integer :: xpe(numFL,numLP,2**ndim)
458  integer :: xgd(numFL,numLP,2**ndim)
459  double precision :: Tcoff_line(numFL)
460  double precision :: weightIndex,weight,ds
461  integer :: i,j,k,igrid,iigrid,ixO^L,ixc^L,ixc^D
462  double precision :: dxMax^D,dxb^D
463 
464  !> interpolate lines into grid cells
465  weightindex=2.d0
466  dxmax^d=2.d0*dl;
467  ixo^l=ixm^ll;
468  do i=1,numfl
469  do j=1,numr(i)
470  do k=1,2**ndim
471  if(mype .eq. xpe(i,j,k)) then
472  igrid=xgd(i,j,k)
473  if(igrid .le. igrids(igridstail)) then
474  ^d&dxb^d=rnode(rpdx^d_,igrid)\
475  ^d&ixcmin^d=floor((xf(i,j,^d)-dxmax^d-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d\
476  ^d&ixcmax^d=floor((xf(i,j,^d)+dxmax^d-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d\
477  {if (ixcmin^d<ixomin^d) ixcmin^d=ixomin^d\}
478  {if (ixcmax^d>ixomax^d) ixcmax^d=ixomax^d\}
479  {do ixc^d=ixcmin^d,ixcmax^d\}
480  ds=0.d0
481  {ds=ds+(xf(i,j,^d)-ps(igrid)%x(ixc^dd,^d))**2\}
482  ds=sqrt(ds)
483  if(ds .le. 0.099d0*dl) then
484  weight=(1/(0.099d0*dl))**weightindex
485  else
486  weight=(1/ds)**weightindex
487  endif
488  ps(igrid)%wextra(ixc^d,tweight_)=ps(igrid)%wextra(ixc^d,tweight_)+weight
489  ps(igrid)%wextra(ixc^d,tcoff_)=ps(igrid)%wextra(ixc^d,tcoff_)+weight*tcoff_line(i)
490  {enddo\}
491  else
492  call mpistop("we need to check here 366Line in mod_trac.t")
493  end if
494  end if
495  end do
496  end do
497  end do
498  ! finish interpolation
499  do iigrid=1,igridstail; igrid=igrids(iigrid);
500  where (ps(igrid)%wextra(ixo^s,tweight_)>0.d0)
501  ps(igrid)%wextra(ixo^s,tcoff_)=ps(igrid)%wextra(ixo^s,tcoff_)/ps(igrid)%wextra(ixo^s,tweight_)
502  elsewhere
503  ps(igrid)%wextra(ixo^s,tcoff_)=0.2d0*tmax
504  endwhere
505  enddo
506  end subroutine block_interp_grid
507 
508  ! TODO remove, not used?
509 ! subroutine get_Tmax_grid(x,w,ixI^L,ixO^L,Tmax_grid)
510 ! integer, intent(in) :: ixI^L,ixO^L
511 ! double precision, intent(in) :: x(ixI^S,1:ndim)
512 ! double precision, intent(out) :: w(ixI^S,1:nw)
513 ! double precision :: Tmax_grid
514 ! double precision :: pth(ixI^S),Te(ixI^S)
515 !
516 ! call mhd_get_pthermal(w,x,ixI^L,ixO^L,pth)
517 ! Te(ixO^S)=pth(ixO^S)/w(ixO^S,rho_)
518 ! Tmax_grid=maxval(Te(ixO^S))
519 ! end subroutine get_Tmax_grid
520 
521  subroutine init_trac_tcoff()
522  integer :: ixI^L,ixO^L,igrid,iigrid
523 
524  ixi^l=ixg^ll;
525  ixo^l=ixm^ll;
526 
527  do iigrid=1,igridstail; igrid=igrids(iigrid);
528  ps(igrid)%wextra(ixi^s,tcoff_)=0.d0
529  ps(igrid)%wextra(ixi^s,tweight_)=0.d0
530  end do
531 
532  end subroutine init_trac_tcoff
533 
534  subroutine update_pegrid()
535 
536  ngrid_trac=0
537  ngrid_ground=0
538  trac_pe(:)=.false.
539  trac_grid(:)=0
540  ground_grid(:)=0
541 
542  ! traverse gridtable to information of grid inside mask region
543  call traverse_gridtable()
544 
545  end subroutine update_pegrid
546 
547  subroutine traverse_gridtable()
549 
550  double precision :: dxb^D,xb^L
551  integer :: iigrid,igrid,j
552  logical, allocatable :: trac_pe_recv(:)
553  double precision :: hcmax_bt
554 
555  allocate(trac_pe_recv(npe))
556 
557  hcmax_bt=xprobmin^nd+(xprobmax^nd-xprobmin^nd)/(domain_nx^nd*2**(refine_max_level-1))
558 
559  do iigrid=1,igridstail; igrid=igrids(iigrid);
560  xbmin^nd=rnode(rpxmin^nd_,igrid)
561  if (xbmin^nd<phys_trac_mask) then
562  trac_pe(mype)=.true.
564  trac_grid(ngrid_trac)=igrid
565  if (xbmin^nd<hcmax_bt) then
568  endif
569  endif
570  enddo
571 
572  call mpi_allreduce(trac_pe,trac_pe_recv,npe,mpi_logical,mpi_lor,icomm,ierrmpi)
573  trac_pe=trac_pe_recv
574 
575  deallocate(trac_pe_recv)
576 
577  end subroutine traverse_gridtable
578 
579  subroutine get_te_grid()
580  integer :: ixI^L,ixO^L,igrid,iigrid,j
581 
582  ixi^l=ixg^ll;
583  ixo^l=ixm^ll;
584 
585  do j=1,ngrid_trac
586  igrid=trac_grid(j)
587 
588  call mhd_get_pthermal(ps(igrid)%w,ps(igrid)%x,ixi^l,ixi^l,ps(igrid)%wextra(ixi^s,tcoff_))
589  !TODO move outside loop for optimziation?
590  if(has_equi_rho0) then
591  ps(igrid)%wextra(ixi^s,tcoff_)=ps(igrid)%wextra(ixi^s,tcoff_)/&
592  (ps(igrid)%w(ixi^s,rho_) + ps(igrid)%equi_vars(ixi^s,equi_rho0_,0))
593  else
594  ps(igrid)%wextra(ixi^s,tcoff_)=ps(igrid)%wextra(ixi^s,tcoff_)/ps(igrid)%w(ixi^s,rho_)
595  endif
596  enddo
597 
598  end subroutine get_te_grid
599 
600  subroutine get_btracing_dir(ipel,igridl,forwardl)
601  integer :: ipel(numFL,numLP),igridl(numFL,numLP)
602  logical :: forwardl(numFL)
603 
604  integer :: igrid,ixO^L,iFL,j,ix^D,idir,ixb^D,ixbb^D
605  double precision :: xb^L,dxb^D,xd^D,factor,Bh
606  integer :: numL(ndim),ixmin(ndim),ixmax(ndim),ix(ndim)
607  logical :: forwardRC(numFL)
608 
609  ipel=-1
610  igridl=0
611  forwardl=.true.
612 
613  ^d&ixomin^d=ixmlo^d;
614  ^d&ixomax^d=ixmhi^d;
615 
616  do j=1,ngrid_ground
617  igrid=ground_grid(j)
618  ^d&dxb^d=rnode(rpdx^d_,igrid);
619  ^d&xbmin^d=rnode(rpxmin^d_,igrid);
620  ^d&xbmax^d=rnode(rpxmax^d_,igrid);
621  ^d&ixmin(^d)=floor((xbmin^d-xprobmin^d)/(mhd_trac_finegrid*dl))+1;
622  ^d&ixmax(^d)=floor((xbmax^d-xprobmin^d)/(mhd_trac_finegrid*dl));
623  ^d&numl(^d)=floor((xprobmax^d-xprobmin^d)/(mhd_trac_finegrid*dl));
624  ixmin(^nd)=1
625  ixmax(^nd)=1
626  numl(^nd)=1
627 
628  {do ix^db=ixmin(^db),ixmax(^db)\}
629  ^d&ix(^d)=ix^d;
630  ifl=0
631  do idir=ndim-1,1,-1
632  ifl=ifl+(ix(idir)-(ndim-1-idir))*(numl(idir))**(ndim-1-idir)
633  enddo
634  ipel(ifl,1)=mype
635  igridl(ifl,1)=igrid
636 
637  ^d&ixb^d=floor((xfi(ifl,^d)-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d;
638  ^d&xd^d=(xfi(ifl,^d)-ps(igrid)%x(ixb^dd,^d))/dxb^d;
639  bh=0.d0
640  {do ixbb^d=0,1\}
641  factor={abs(1-ix^d-xd^d)*}
642  if (b0field) then
643  bh=bh+factor*(ps(igrid)%w(ixb^d+ixbb^d,mag(^nd))+ps(igrid)%B0(ixb^d+ixbb^d,^nd,0))
644  else
645  bh=bh+factor*ps(igrid)%w(ixb^d+ixbb^d,mag(^nd))
646  endif
647  {enddo\}
648  if (bh>0) then
649  forwardl(ifl)=.true.
650  else
651  forwardl(ifl)=.false.
652  endif
653  {enddo\}
654  enddo
655 
656  call mpi_allreduce(forwardl,forwardrc,numfl,mpi_logical,&
657  mpi_land,icomm,ierrmpi)
658  forwardl=forwardrc
659 
660  end subroutine get_btracing_dir
661 
662  subroutine get_tcoff_line(xFL,numR,TcoffFL,ipeFL,igridFL,forwardFL,mask)
663  use mod_trace_field
664 
665  double precision :: xFL(numFL,numLP,ndim)
666  integer :: numR(numFL)
667  double precision :: TcoffFL(numFL),TmaxFL(numFL)
668  integer :: ipeFL(numFL,numLP),igridFL(numFL,numLP)
669  logical :: forwardFL(numFL)
670  logical, intent(in) :: mask
671 
672  integer :: nwP,nwL,iFL,iLP
673  double precision :: wPm(numFL,numLP,2),wLm(numFL,1+2)
674  character(len=std_len) :: ftype,tcondi
675 
676  nwp=2
677  nwl=2
678  ftype='Bfield'
679  tcondi='TRAC'
680  call trace_field_multi(xfl,wpm,wlm,dl,numfl,numlp,nwp,nwl,forwardfl,ftype,tcondi)
681  do ifl=1,numfl
682  numr(ifl)=int(wlm(ifl,1))
683  tcofffl(ifl)=wlm(ifl,2)
684  tmaxfl(ifl)=wlm(ifl,3)
685  if(mask) then
686  if(tcofffl(ifl)>0.2d0*tmax) tcofffl(ifl)=0.2d0*tmax
687  else
688  tmaxfl(ifl)=wlm(ifl,3)
689  if(tcofffl(ifl)>0.2d0*tmaxfl(ifl)) tcofffl(ifl)=0.2d0*tmaxfl(ifl)
690  end if
691 
692  if(tcofffl(ifl)<t_bott) tcofffl(ifl)=t_bott
693  enddo
694 
695  do ifl=1,numfl
696  if (numr(ifl)>0) then
697  do ilp=1,numr(ifl)
698  ipefl(ifl,ilp)=int(wpm(ifl,ilp,1))
699  igridfl(ifl,ilp)=int(wpm(ifl,ilp,2))
700  enddo
701  endif
702  enddo
703 
704 
705  end subroutine get_tcoff_line
706 
707  subroutine interp_tcoff(xF,ipel,igridl,numR,Tlcoff)
708  double precision :: xF(numFL,numLP,ndim)
709  integer :: numR(numFL),ipel(numFL,numLP),igridl(numFL,numLP)
710  double precision :: Tlcoff(numFL)
711 
712  integer :: iFL,iLP,ixO^L,ixI^L,ixc^L,ixb^L,ixc^D
713  integer :: igrid,j,ipmin,ipmax,igrid_nb
714  double precision :: dxb^D,dxMax^D,xb^L,Tcnow
715  double precision :: xFnow(ndim)
716  integer :: weightIndex,idn^D,ixmax^ND
717  double precision :: ds,weight
718 
719  weightindex=2
720 
721  ^d&iximin^d=ixglo^d;
722  ^d&iximax^d=ixghi^d;
723  ^d&ixomin^d=ixmlo^d;
724  ^d&ixomax^d=ixmhi^d;
725 
726  do ifl=1,numfl
727  ilp=1
728  tcnow=tlcoff(ifl)
729  do while (ilp<=numr(ifl))
730  ! find points in the same grid
731  ipmin=ilp
732  do while (ipel(ifl,ipmin)/=mype .and. ipmin<=numr(ifl))
733  ipmin=ipmin+1
734  enddo
735  igrid=igridl(ifl,ipmin)
736  ipmax=ipmin
737  do while (ipel(ifl,ipmax)==mype .and. igridl(ifl,ipmax+1)==igrid .and. ipmax<numr(ifl))
738  ipmax=ipmax+1
739  enddo
740 
741  ! parameters for the grid
742  ^d&dxb^d=rnode(rpdx^d_,igrid);
743  ^d&dxmax^d=4*dxb^d;
744 
745  do ilp=ipmin,ipmax
746  xfnow(:)=xf(ifl,ilp,:)
747  ^d&ixbmin^d=floor((xfnow(^d)-dxmax^d-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d;
748  ^d&ixbmax^d=floor((xfnow(^d)+dxmax^d-ps(igrid)%x(ixomin^dd,^d))/dxb^d)+ixomin^d;
749 
750  ! do interpolation inside the grid
751  {ixcmin^d=max(ixbmin^d,ixomin^d)\}
752  {ixcmax^d=min(ixbmax^d,ixomax^d)\}
753  xbmin^nd=rnode(rpxmin^nd_,igrid)
754  xbmax^nd=rnode(rpxmax^nd_,igrid)
755  ixmax^nd=floor((phys_trac_mask-xbmin^nd)/dxb^nd)+ixomin^nd
756  if (xbmax^nd>phys_trac_mask) ixcmax^nd=min(ixmax^nd,ixcmax^nd)
757  {do ixc^d=ixcmin^d,ixcmax^d\}
758  ds=0.d0
759  {ds=ds+(xfnow(^d)-ps(igrid)%x(ixc^dd,^d))**2\}
760  ds=sqrt(ds)
761  if(ds<1.0d-2*dxb1) then
762  weight=(1/(1.0d-2*dxb1))**weightindex
763  else
764  weight=(1/ds)**weightindex
765  endif
766  ps(igrid)%wextra(ixc^d,tweight_)=ps(igrid)%wextra(ixc^d,tweight_)+weight
767  ps(igrid)%wextra(ixc^d,tcoff_)=ps(igrid)%wextra(ixc^d,tcoff_)+weight*tcnow
768  {enddo\}
769 
770  ! do interpolation in neighbor grids that have the same level and pe
771  {
772  if (ixbmin^d<ixomin^d) then
773  idn^dd=0;
774  idn^d=-1
775  if (neighbor(2,idn^dd,igrid)==mype .and. neighbor_type(idn^dd,igrid)==neighbor_sibling) then
776  igrid_nb=neighbor(1,idn^dd,igrid)
777  ixcmin^dd=max(ixbmin^dd,ixomin^dd);
778  ixcmax^dd=min(ixbmax^dd,ixomax^dd);
779  ixcmin^d=ixomax^d+(ixbmin^d-ixomin^d)
780  ixcmax^d=ixomax^d
781  xbmin^nd=rnode(rpxmin^nd_,igrid_nb)
782  xbmax^nd=rnode(rpxmax^nd_,igrid_nb)
783  ixmax^nd=floor((phys_trac_mask-xbmin^nd)/dxb^nd)+ixomin^nd
784  if (xbmax^nd>phys_trac_mask) ixcmax^nd=min(ixmax^nd,ixcmax^nd)
785 
786  {do ixc^dd=ixcmin^dd,ixcmax^dd;}
787  ds=0.d0
788  {ds=ds+(xfnow(^dd)-ps(igrid_nb)%x({ixc^dd},^dd))**2;}
789  ds=sqrt(ds)
790  if(ds<1.0d-2*dxb1) then
791  weight=(1/(1.0d-2*dxb1))**weightindex
792  else
793  weight=(1/ds)**weightindex
794  endif
795  ps(igrid_nb)%wextra(ixc^dd,tweight_)=ps(igrid_nb)%wextra(ixc^dd,tweight_)+weight
796  ps(igrid_nb)%wextra(ixc^dd,tcoff_)=ps(igrid_nb)%wextra(ixc^dd,tcoff_)+weight*tcnow
797  {enddo;}
798  endif
799  endif
800 
801  if (ixbmax^d>ixomin^d) then
802  idn^dd=0;
803  idn^d=1
804  if (neighbor(2,idn^dd,igrid)==mype .and. neighbor_type(idn^dd,igrid)==neighbor_sibling) then
805  igrid_nb=neighbor(1,idn^dd,igrid)
806  xbmin^nd=rnode(rpxmin^nd_,igrid_nb)
807  if (xbmin^nd<phys_trac_mask) then
808  ixcmin^dd=max(ixbmin^dd,ixomin^dd);
809  ixcmax^dd=min(ixbmax^dd,ixomax^dd);
810  ixcmin^d=ixomin^d
811  ixcmax^d=ixomin^d+(ixbmax^d-ixomax^d)
812  xbmax^nd=rnode(rpxmax^nd_,igrid_nb)
813  ixmax^nd=floor((phys_trac_mask-xbmin^nd)/dxb^nd)+ixomin^nd
814  if (xbmax^nd>phys_trac_mask) ixcmax^nd=min(ixmax^nd,ixcmax^nd)
815 
816  {do ixc^dd=ixcmin^dd,ixcmax^dd;}
817  ds=0.d0
818  {ds=ds+(xfnow(^dd)-ps(igrid_nb)%x({ixc^dd},^dd))**2;}
819  ds=sqrt(ds)
820  if(ds<1.0d-2*dxb1) then
821  weight=(1/(1.0d-2*dxb1))**weightindex
822  else
823  weight=(1/ds)**weightindex
824  endif
825  ps(igrid_nb)%wextra(ixc^dd,tweight_)=ps(igrid_nb)%wextra(ixc^dd,tweight_)+weight
826  ps(igrid_nb)%wextra(ixc^dd,tcoff_)=ps(igrid_nb)%wextra(ixc^dd,tcoff_)+weight*tcnow
827  {enddo;}
828  endif
829  endif
830  endif
831  \}
832  enddo
833 
834  ilp=ipmax+1
835  enddo
836  enddo
837 
838 
839  do j=1,ngrid_trac
840  igrid=trac_grid(j)
841  where(ps(igrid)%wextra(ixo^s,tweight_)>0.d0)
842  ps(igrid)%wextra(ixo^s,tcoff_)=ps(igrid)%wextra(ixo^s,tcoff_)/ps(igrid)%wextra(ixo^s,tweight_)
843  elsewhere
844  ps(igrid)%wextra(ixo^s,tcoff_)=t_bott
845  endwhere
846  enddo
847 
848  end subroutine interp_tcoff
849 
850 end module mod_trac
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 domain_nx
number of cells for each dimension in level-one mesh
double precision phys_trac_mask
integer, parameter ndim
Number of spatial dimensions for grid variables.
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.
integer refine_max_level
Maximal number of AMR levels.
integer max_blocks
The maximum number of grid blocks in a processor.
Definition: mod_mhd.t:1
double precision dxt
Definition: mod_trac.t:11
subroutine get_btracing_dir(ipel, igridl, forwardl)
Definition: mod_trac.t:601
double precision t_bott
Definition: mod_trac.t:7
integer ngrid_trac
Definition: mod_trac.t:16
subroutine get_te_grid()
Definition: mod_trac.t:580
integer numlp
Definition: mod_trac.t:6
subroutine rk_bdir(nowgridc, nowpoint, ix_nextD, now_dir, Bdir, ix_mod, first, init_dir)
Definition: mod_trac.t:405
integer ngrid_ground
Definition: mod_trac.t:16
subroutine tracb(mask, T_peak)
Definition: mod_trac.t:219
subroutine init_trac_block(mask)
Definition: mod_trac.t:78
double precision, dimension(:,:), allocatable xfi
Definition: mod_trac.t:8
integer numfl
Definition: mod_trac.t:6
double precision, dimension(:^d &,:), allocatable xt
Definition: mod_trac.t:13
double precision, dimension(ndim) xtmax
Definition: mod_trac.t:12
subroutine interp_tcoff(xF, ipel, igridl, numR, Tlcoff)
Definition: mod_trac.t:708
integer d
Definition: mod_trac.t:10
subroutine block_estable(mask, Tcoff, Tcmax, Bdir, peArr, gdArr)
Definition: mod_trac.t:242
integer numxt
Definition: mod_trac.t:10
subroutine update_pegrid()
Definition: mod_trac.t:535
subroutine get_tcoff_line(xFL, numR, TcoffFL, ipeFL, igridFL, forwardFL, mask)
Definition: mod_trac.t:663
subroutine block_trace_mfl(mask, Tcoff, Tcoff_line, Tcmax, Bdir, peArr, gdArr, xF, numR, xpe, xgd)
Definition: mod_trac.t:303
logical, dimension(:), allocatable trac_pe
Definition: mod_trac.t:17
integer, dimension(:), allocatable ground_grid
Definition: mod_trac.t:15
double precision dl
Definition: mod_trac.t:7
double precision trac_delta
Definition: mod_trac.t:7
subroutine ltrac(T_peak)
Definition: mod_trac.t:165
double precision tmax
Definition: mod_trac.t:7
subroutine block_interp_grid(mask, xF, numR, xpe, xgd, Tcoff_line)
Definition: mod_trac.t:454
subroutine init_trac_tcoff()
Definition: mod_trac.t:522
integer, dimension(:), allocatable trac_grid
Definition: mod_trac.t:15
subroutine traverse_gridtable()
Definition: mod_trac.t:548
subroutine trac_simple(tco_global, trac_alfa, T_peak)
Definition: mod_trac.t:143
subroutine tracl(mask, T_peak)
Definition: mod_trac.t:181
double precision, dimension(ndim) xtmin
Definition: mod_trac.t:12
subroutine init_trac_line(mask)
Definition: mod_trac.t:21
subroutine trace_field_multi(xfm, wPm, wLm, dL, numL, numP, nwP, nwL, forwardm, ftype, tcondi)