MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_trac.t
Go to the documentation of this file.
1module mod_trac
4 use mod_comm_lib, only: mpistop
5 implicit none
6 private
7 ! common
8 integer :: numFL,numLP
9 double precision :: dL,Tmax,trac_delta,T_bott
10 double precision, allocatable :: xFi(:,:)
11 ! TRACB
12 integer :: numxT^D
13 double precision :: dxT^D
14 double precision :: xTmin(ndim),xTmax(ndim)
15 double precision, allocatable :: xT(:^D&,:)
16 ! TRACL mask
17 integer, allocatable :: trac_grid(:),ground_grid(:)
18 integer :: ngrid_trac,ngrid_ground
19 logical, allocatable :: trac_pe(:)
21contains
22
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)
28 integer :: memxFi
29
30 trac_delta=0.25d0
31 !----------------- init magnetic field -------------------!
32 refine_factor=2**(refine_max_level-1)
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)
37 ! max length of a field line
38 if (.not. mask) phys_trac_mask=xprobmax^nd
39 {^iftwod
40 lengthfl=(phys_trac_mask-xprobmin(2))*3.d0
41 }
42 {^ifthreed
43 lengthfl=(phys_trac_mask-xprobmin(3))*3.d0
44 }
45
46 numlp=floor(lengthfl/dl)
47 numl=1
48 numfl=1
49 do j=1,ndim-1
50 ! number of field lines, every 4 finest grid, every direction
51 finegrid=phys_trac_finegrid
52 numl(j)=floor((xprobmax(j)-xprobmin(j))/dl/finegrid)
53 numfl=numfl*numl(j)
54 end do
55 allocate(xfi(numfl,ndim))
56 xfi(:,ndim)=xprobmin(ndim)+dl/50.d0
57 {do ix^db=1,numl(^db)\}
58 ^d&ix(^d)=ix^d\
59 ifl=0
60 do j=ndim-1,1,-1
61 ifl=ifl+(ix(j)-(ndim-1-j))*(numl(j))**(ndim-1-j)
62 end do
63 xfi(ifl,1:ndim-1)=xprobmin(1:ndim-1)+finegrid*dl*ix(1:ndim-1)-finegrid*dl/2.d0
64 {end do\}
65 {^iftwod
66 if(mype .eq. 0) write(*,*) 'NOTE: 2D TRAC method take the y-dir == grav-dir'
67 }
68 {^ifthreed
69 if(mype .eq. 0) write(*,*) 'NOTE: 3D TRAC method take the z-dir == grav-dir'
70 }
71
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'
75
76 allocate(trac_pe(0:npe-1),trac_grid(max_blocks),ground_grid(max_blocks))
77
78 end subroutine init_trac_line
79
80 subroutine init_trac_block(mask)
81 logical, intent(in) :: mask
82 integer :: refine_factor,finegrid,iFL,j
83 integer :: ix^D,ixT^L
84 integer :: numL(ndim),ix(ndim)
85 double precision :: lengthFL
86 double precision :: ration,a0
87 double precision :: xprobmin(ndim),xprobmax(ndim),dxT(ndim)
88
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)\
93 ^d&dxt(ndim)=dxt^d\
94 finegrid=phys_trac_finegrid
95 {^ifoned
96 dl=dxt1/finegrid
97 }
98 {^nooned
99 dl=min(dxt^d)/finegrid
100 }
101 ! table for interpolation
102 ^d&xtmin(^d)=xprobmin^d\
103 ^d&xtmax(^d)=xprobmax^d\
104 if(mask) xtmax(ndim)=phys_trac_mask
105 ! max length of a field line
106 if(mask) then
107 lengthfl=maxval(xtmax-xprobmin)*3.d0
108 else
109 lengthfl=maxval(xprobmax-xprobmin)*3.d0
110 end if
111 numlp=floor(lengthfl/dl)
112 ^d&numxt^d=ceiling((xtmax(^d)-xtmin(^d)-smalldouble)/dxt^d)\
113 allocate(xt(numxt^d,ndim))
114 ixtmin^d=1\
115 ixtmax^d=numxt^d\
116 {do j=1,numxt^d
117 xt(j^d%ixT^s,^d)=(j-0.5d0)*dxt^d+xtmin(^d)
118 end do\}
119 if(mask) xtmax(ndim)=maxval(xt(:^d&,ndim))+half*dxt(ndim)
120 numl=1
121 numfl=1
122 do j=1,ndim-1
123 ! number of field lines, every 4 finest grid, every direction
124 numl(j)=floor((xprobmax(j)-xprobmin(j))/dl)
125 numfl=numfl*numl(j)
126 end do
127 allocate(xfi(numfl,ndim))
128 xfi(:,ndim)=xprobmin(ndim)+dl/50.d0
129 {do ix^db=1,numl(^db)\}
130 ^d&ix(^d)=ix^d\
131 ifl=0
132 do j=ndim-1,1,-1
133 ifl=ifl+(ix(j)-(ndim-1-j))*(numl(j))**(ndim-1-j)
134 end do
135 xfi(ifl,1:ndim-1)=xprobmin(1:ndim-1)+dl*ix(1:ndim-1)-dl/2.d0
136 {end do\}
137 {^iftwod
138 if(mype .eq. 0) write(*,*) 'NOTE: 2D TRAC method take the y-dir == grav-dir'
139 }
140 {^ifthreed
141 if(mype .eq. 0) write(*,*) 'NOTE: 3D TRAC method take the z-dir == grav-dir'
142 }
143 end subroutine init_trac_block
144
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
148
149 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
150 {^ifoned
151 ps(igrid)%special_values(1)=tco_global
152 }
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)
155 end if
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
160 end if
161 ps(igrid)%wextra(ixg^t,iw_tcoff)=ps(igrid)%special_values(1)
162 !> special values(2) to save old tcutoff
163 ps(igrid)%special_values(2)=ps(igrid)%special_values(1)
164 end do
165 end subroutine trac_simple
166
167 subroutine ltrac(T_peak)
168 double precision, intent(in) :: T_peak
169 integer :: iigrid, igrid
170 integer :: ixO^L,trac_tcoff
171
172 ixo^l=ixm^ll^ladd1;
173 trac_tcoff=iw_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
179 end where
180 end do
181 end subroutine ltrac
182
183 subroutine tracl(mask,T_peak)
184 logical, intent(in) :: mask
185 double precision, intent(in) :: T_peak
186 integer :: ix^D,j
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)
193
194 tmax=t_peak
195 xf=zero
196
197 do j=1,ndim
198 xf(:,1,j)=xfi(:,j)
199 end do
200
201 call mpi_barrier(icomm,ierrmpi)
202 ! record pe and grid for mask region
203 call update_pegrid()
204 ! get temperature for the calculation of Te gradient,
205 ! which will be stored in wextra(ixI^S,iw_Tcoff) temporarily
206 call get_te_grid()
207 ! find out direction for tracing B field
208 call get_btracing_dir(ipel,igridl,forwardl)
209 ! trace Bfield and get Tcoff for each field line
210 call get_tcoff_line(xf,numr,tlcoff,ipel,igridl,forwardl,mask)
211 ! init vairable Tcoff_ and Tweight_
212 call init_trac_tcoff()
213 ! get cell Tcoff via interpolation
214 call interp_tcoff(xf,ipel,igridl,numr,tlcoff)
215 ! convert primitive data back to conserved data
216 call mpi_barrier(icomm,ierrmpi)
217
218
219 end subroutine tracl
220
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)
229
230 tmax=t_peak
231 tcoff=zero
232 tcmax=zero
233 bdir=zero
234 pearr=-1
235 gdarr=-1
236 call block_estable(mask,tcoff,tcmax,bdir,pearr,gdarr)
237 xf=zero
238 numr=0
239 tcoff_line=zero
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)
242 end subroutine tracb
243
244 subroutine block_estable(mask,Tcoff,Tcmax,Bdir,peArr,gdArr)
245 logical :: mask
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
254
255 tcoff_recv=zero
256 tcmax_recv=zero
257 bdir_recv=zero
258 pearr_recv=-1
259 gdarr_recv=-1
260 !> combine table from different processors
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 \}
271 intab=0
272 {if (ix^d .le. numxt^d) intab=intab+1 \}
273 if(intab .eq. ndim) then
274 !> in principle, no overlap will happen here
275 tcoff(ix^d)=max(tcoff(ix^d),ps(igrid)%special_values(1))
276 tcmax(ix^d)=ps(igrid)%special_values(2)
277 !> abs(Bdir) <= 1, so that Bdir+2 should always be positive
278 bdir(ix^d,1:ndim)=ps(igrid)%special_values(3:3+ndim-1)+2.d0
279 pearr(ix^d)=mype
280 gdarr(ix^d)=igrid
281 end if
282 {end do\}
283 end do
284 call mpi_barrier(icomm,ierrmpi)
285 numxt={numxt^d*}
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)
294 if(.not. mask) then
295 call mpi_allreduce(tcmax,tcmax_recv,numxt,mpi_double_precision,&
296 mpi_max,icomm,ierrmpi)
297 end if
298 pearr=pearr_recv
299 gdarr=gdarr_recv
300 tcoff=tcoff_recv
301 bdir=bdir_recv-2.d0
302 if(.not. mask) tcmax=tcmax_recv
303 end subroutine block_estable
304
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)
318
319 do i=1,numfl
320 flag=.true.
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)
324 ix_next^d=k^d\
325 j=1
326 xf(i,j,:)=xfi(i,:)
327 do while(flag)
328 nowpoint(:)=xf(i,j,:)
329 nowgridc(:)=xt(ix_next^d,:)
330 first=.true.
331 if(j .eq. 1) then
332 call rk_bdir(nowgridc,nowpoint,ix_next^d,now_dir1,bdir,&
333 ix_mod,first,init_dir)
334 else
335 call rk_bdir(nowgridc,nowpoint,ix_next^d,now_dir1,bdir,&
336 ix_mod,first)
337 end if
338 {^iftwod
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))
347 }
348 {^ifthreed
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))
365 }
366 nowpoint(:)=nowpoint(:)+init_dir*now_dir1*dl
367 {if(nowpoint(^d) .gt. xtmax(^d) .or. nowpoint(^d) .lt. xtmin(^d)) then
368 flag=.false.
369 end if\}
370 if(mask .and. nowpoint(ndim) .gt. phys_trac_mask) then
371 flag=.false.
372 end if
373 if(flag) then
374 first=.false.
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,&
378 ix_mod,first)
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
381 flag=.false.
382 end if\}
383 if(mask .and. xf(i,j+1,ndim) .gt. phys_trac_mask) then
384 flag=.false.
385 end if
386 if(flag) then
387 ^d&ix_next^d=ceiling((xf(i,j+1,^d)-xtmin(^d)-smalldouble)/dxt^d)\
388 j=j+1
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))
391 end if
392 end if
393 end do
394 numr(i)=j
395 if(mask) then
396 if(tcoff_line(i) .gt. tmax*0.2d0) then
397 tcoff_line(i)=tmax*0.2d0
398 end if
399 else
400 if(tcoff_line(i) .gt. tcmax_line(i)*0.2d0) then
401 tcoff_line(i)=tcmax_line(i)*0.2d0
402 end if
403 end if
404 end do
405 end subroutine block_trace_mfl
406
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)
414 logical :: first
415 double precision, optional :: init_dir
416
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)
423 else
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)
428 end if
429 else
430 ix_mod(^d,:)=ix_next^d
431 alfa_mod(^d,:)=half
432 end if\}
433 now_dir=zero
434 {^iftwod
435 do k1=1,2
436 do k2=1,2
437 now_dir=now_dir + bdir(ix_mod(1,k1),ix_mod(2,k2),:)*alfa_mod(1,k1)*alfa_mod(2,k2)
438 end do
439 end do
440 }
441 {^ifthreed
442 do k1=1,2
443 do k2=1,2
444 do k3=1,2
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)
447 end do
448 end do
449 end do
450 }
451 if(present(init_dir)) then
452 init_dir=sign(one,now_dir(ndim))
453 end if
454 end subroutine rk_bdir
455
456 subroutine block_interp_grid(mask,xF,numR,xpe,xgd,Tcoff_line)
457 logical :: mask
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
466
467 !> interpolate lines into grid cells
468 weightindex=2.d0
469 dxmax^d=2.d0*dl;
470 ixo^l=ixm^ll;
471 do i=1,numfl
472 do j=1,numr(i)
473 do k=1,2**ndim
474 if(mype .eq. xpe(i,j,k)) then
475 igrid=xgd(i,j,k)
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\}
483 ds=0.d0
484 {ds=ds+(xf(i,j,^d)-ps(igrid)%x(ixc^dd,^d))**2\}
485 ds=sqrt(ds)
486 if(ds .le. 0.099d0*dl) then
487 weight=(1/(0.099d0*dl))**weightindex
488 else
489 weight=(1/ds)**weightindex
490 endif
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)
493 {enddo\}
494 else
495 call mpistop("we need to check here 366Line in mod_trac.t")
496 end if
497 end if
498 end do
499 end do
500 end do
501 ! finish interpolation
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)
505 elsewhere
506 ps(igrid)%wextra(ixo^s,iw_tcoff)=0.2d0*tmax
507 endwhere
508 enddo
509 end subroutine block_interp_grid
510
511 ! TODO remove, not used?
512! subroutine get_Tmax_grid(x,w,ixI^L,ixO^L,Tmax_grid)
513! integer, intent(in) :: ixI^L,ixO^L
514! double precision, intent(in) :: x(ixI^S,1:ndim)
515! double precision, intent(out) :: w(ixI^S,1:nw)
516! double precision :: Tmax_grid
517! double precision :: pth(ixI^S),Te(ixI^S)
518!
519! call mhd_get_pthermal(w,x,ixI^L,ixO^L,pth)
520! Te(ixO^S)=pth(ixO^S)/w(ixO^S,rho_)
521! Tmax_grid=maxval(Te(ixO^S))
522! end subroutine get_Tmax_grid
523
524 subroutine init_trac_tcoff()
525 integer :: ixI^L,ixO^L,igrid,iigrid
526
527 ixi^l=ixg^ll;
528 ixo^l=ixm^ll;
529
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
533 end do
534 end subroutine init_trac_tcoff
535
536 subroutine update_pegrid()
537
538 ngrid_trac=0
539 ngrid_ground=0
540 trac_pe(:)=.false.
541 trac_grid(:)=0
542 ground_grid(:)=0
543
544 ! traverse gridtable to information of grid inside mask region
545 call traverse_gridtable()
546
547 end subroutine update_pegrid
548
549 subroutine traverse_gridtable()
551
552 double precision :: dxb^D,xb^L
553 integer :: iigrid,igrid,j
554 logical, allocatable :: trac_pe_recv(:)
555 double precision :: hcmax_bt
556
557 allocate(trac_pe_recv(npe))
558
559 hcmax_bt=xprobmin^nd+(xprobmax^nd-xprobmin^nd)/(domain_nx^nd*2**(refine_max_level-1))
560
561 do iigrid=1,igridstail; igrid=igrids(iigrid);
562 xbmin^nd=rnode(rpxmin^nd_,igrid)
563 if (xbmin^nd<phys_trac_mask) then
564 trac_pe(mype)=.true.
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
570 endif
571 endif
572 enddo
573
574 call mpi_allreduce(trac_pe,trac_pe_recv,npe,mpi_logical,mpi_lor,icomm,ierrmpi)
575 trac_pe=trac_pe_recv
576
577 deallocate(trac_pe_recv)
578 end subroutine traverse_gridtable
579
580 subroutine get_te_grid()
581 integer :: ixI^L,ixO^L,igrid,iigrid,j
582 double precision :: rho(ixM^T)
583
584 ixi^l=ixg^ll;
585 ixo^l=ixm^ll;
586 do j=1,ngrid_trac
587 igrid=trac_grid(j)
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)
591 enddo
592 end subroutine get_te_grid
593
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)
601
602 ipel=-1
603 igridl=0
604 forwardl=.true.
605 ^d&ixomin^d=ixmlo^d;
606 ^d&ixomax^d=ixmhi^d;
607
608 do j=1,ngrid_ground
609 igrid=ground_grid(j)
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));
616 ixmin(^nd)=1
617 ixmax(^nd)=1
618 numl(^nd)=1
619 {do ix^db=ixmin(^db),ixmax(^db)\}
620 ^d&ix(^d)=ix^d;
621 ifl=0
622 do idir=ndim-1,1,-1
623 ifl=ifl+(ix(idir)-(ndim-1-idir))*(numl(idir))**(ndim-1-idir)
624 enddo
625 ipel(ifl,1)=mype
626 igridl(ifl,1)=igrid
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;
629 bh=0.d0
630 {do ixbb^d=0,1\}
631 factor={abs(1-ix^d-xd^d)*}
632 if (b0field) then
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))
635 else
636 bh=bh+factor*(ps(igrid)%B0(ixb^d+ixbb^d,^nd,0))
637 endif
638 else
639 bh=bh+factor*ps(igrid)%w(ixb^d+ixbb^d,iw_mag(^nd))
640 endif
641 {enddo\}
642 if (bh>0) then
643 forwardl(ifl)=.true.
644 else
645 forwardl(ifl)=.false.
646 endif
647 {enddo\}
648 enddo
649 call mpi_allreduce(forwardl,forwardrc,numfl,mpi_logical,&
650 mpi_land,icomm,ierrmpi)
651 forwardl=forwardrc
652 end subroutine get_btracing_dir
653
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
665
666 nwp=2
667 nwl=2
668 ftype='Bfield'
669 tcondi='TRAC'
670 call trace_field_multi(xfl,wpm,wlm,dl,numfl,numlp,nwp,nwl,forwardfl,ftype,tcondi)
671 do ifl=1,numfl
672 numr(ifl)=int(wlm(ifl,1))
673 tcofffl(ifl)=wlm(ifl,2)
674 tmaxfl(ifl)=wlm(ifl,3)
675 if(mask) then
676 if(tcofffl(ifl)>0.2d0*tmax) tcofffl(ifl)=0.2d0*tmax
677 else
678 tmaxfl(ifl)=wlm(ifl,3)
679 if(tcofffl(ifl)>0.2d0*tmaxfl(ifl)) tcofffl(ifl)=0.2d0*tmaxfl(ifl)
680 end if
681
682 if(tcofffl(ifl)<t_bott) tcofffl(ifl)=t_bott
683 enddo
684
685 do ifl=1,numfl
686 if (numr(ifl)>0) then
687 do ilp=1,numr(ifl)
688 ipefl(ifl,ilp)=int(wpm(ifl,ilp,1))
689 igridfl(ifl,ilp)=int(wpm(ifl,ilp,2))
690 enddo
691 endif
692 enddo
693 end subroutine get_tcoff_line
694
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
705
706 weightindex=2
707
708 ^d&iximin^d=ixglo^d;
709 ^d&iximax^d=ixghi^d;
710 ^d&ixomin^d=ixmlo^d;
711 ^d&ixomax^d=ixmhi^d;
712
713 do ifl=1,numfl
714 ilp=1
715 tcnow=tlcoff(ifl)
716 do while (ilp<=numr(ifl))
717 ! find points in the same grid
718 ipmin=ilp
719 do while (ipel(ifl,ipmin)/=mype .and. ipmin<=numr(ifl))
720 ipmin=ipmin+1
721 enddo
722 igrid=igridl(ifl,ipmin)
723 ipmax=ipmin
724 do while (ipel(ifl,ipmax)==mype .and. igridl(ifl,ipmax+1)==igrid .and. ipmax<numr(ifl))
725 ipmax=ipmax+1
726 enddo
727
728 ! parameters for the grid
729 ^d&dxb^d=rnode(rpdx^d_,igrid);
730 ^d&dxmax^d=4*dxb^d;
731
732 do ilp=ipmin,ipmax
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;
736
737 ! do interpolation inside the grid
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\}
745 ds=0.d0
746 {ds=ds+(xfnow(^d)-ps(igrid)%x(ixc^dd,^d))**2\}
747 ds=sqrt(ds)
748 if(ds<1.0d-2*dxb1) then
749 weight=(1/(1.0d-2*dxb1))**weightindex
750 else
751 weight=(1/ds)**weightindex
752 endif
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
755 {enddo\}
756
757 ! do interpolation in neighbor grids that have the same level and pe
758 {
759 if (ixbmin^d<ixomin^d) then
760 idn^dd=0;
761 idn^d=-1
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)
767 ixcmax^d=ixomax^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)
772
773 {do ixc^dd=ixcmin^dd,ixcmax^dd;}
774 ds=0.d0
775 xc(:)=ps(igrid_nb)%x(ixc^dd,:)
776 {ds=ds+(xfnow(^dd)-xc(^dd))**2;}
777 ds=sqrt(ds)
778 if(ds<1.0d-2*dxb1) then
779 weight=(1/(1.0d-2*dxb1))**weightindex
780 else
781 weight=(1/ds)**weightindex
782 endif
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
785 {enddo;}
786 endif
787 endif
788
789 if (ixbmax^d>ixomin^d) then
790 idn^dd=0;
791 idn^d=1
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);
798 ixcmin^d=ixomin^d
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)
803
804 {do ixc^dd=ixcmin^dd,ixcmax^dd;}
805 ds=0.d0
806 xc(:)=ps(igrid_nb)%x(ixc^dd,:)
807 {ds=ds+(xfnow(^dd)-xc(^dd))**2;}
808 ds=sqrt(ds)
809 if(ds<1.0d-2*dxb1) then
810 weight=(1/(1.0d-2*dxb1))**weightindex
811 else
812 weight=(1/ds)**weightindex
813 endif
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
816 {enddo;}
817 endif
818 endif
819 endif
820 \}
821 enddo
822
823 ilp=ipmax+1
824 enddo
825 enddo
826
827 do j=1,ngrid_trac
828 igrid=trac_grid(j)
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)
831 elsewhere
832 ps(igrid)%wextra(ixo^s,iw_tcoff)=t_bott
833 endwhere
834 enddo
835 end subroutine interp_tcoff
836
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
839 ! default lower limit of cutoff temperature
840 t_bott = t_bott_in
841 select case(phys_trac_type)
842 case(0)
843 !> Test case, fixed cutoff temperature
844 !> do nothing here
845 case(1)
846 !> 1D TRAC method
847 call trac_simple(tco,trac_alfa,t_peak)
848 case(2)
849 !> LTRAC method, by iijima et al. 2021
850 !> do nothing here
851 call ltrac(t_peak)
852 case(3)
853 !> 2D or 3D TRACL(ine) method
854 call tracl(.false.,t_peak)
855 case(4)
856 !> 2D or 3D TRACB(lock) method
857 call tracb(.false.,t_peak)
858 case(5)
859 !> 2D or 3D TRACL(ine) method with mask
860 call tracl(.true.,t_peak)
861 case(6)
862 !> 2D or 3D TRACB(lock) method with mask
863 call tracb(.true.,t_peak)
864 case default
865 call mpistop("undefined TRAC method type")
866 end select
867 end subroutine trac_after_setdt
868
872
873 if(phys_trac) then
874 phys_trac_after_setdt => trac_after_setdt
875 if(phys_trac_type .eq. 3) then
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.)
879 end if
880 if(phys_trac_type .eq. 4) then
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.)
885 end if
886 if(phys_trac_type .eq. 5) then
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.)
890 end if
891 if(phys_trac_type .eq. 6) then
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.)
896 end if
897 end if
898 end subroutine initialize_trac_after_settree
899end module mod_trac
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
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.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
procedure(sub_trac_after_setdt), pointer phys_trac_after_setdt
Definition mod_physics.t:61
subroutine, public initialize_trac_after_settree
Definition mod_trac.t:870
subroutine trace_field_multi(xfm, wpm, wlm, dl, numl, nump, nwp, nwl, forwardm, ftype, tcondi)