19 subroutine getlqgrid(lQgrid,ixI^L,ixO^L,qt,trelax,x,mode,vtim,vlim)
21 integer,
intent(in) :: ixI^L, ixO^L, mode, vtim, vlim
22 double precision,
intent(in) :: qt, trelax, x(ixI^S,1:ndim)
23 double precision,
intent(inout) :: lQgrid(ixI^S)
25 double precision :: tl1,tl2,tl3,tl4,tl5,tl6,t2,stt
26 double precision :: tr1,tr2,tr3,tr4,tr5,tr6
27 double precision :: va0
28 integer :: Bp,Bt1,Bt2,ix^D,i
32 stt = (xprobmax1- xprobmin1) / (dble(vlim) - 1.d0)
37 if(qt .lt. trelax)
then
41 if(qt .lt. trelax +
tb1(i))
then
42 tl1 = dexp(-(t2 -
tb1(i-3))**2/(0.5d0*
ta1(i-3))**2)
43 tl2 = dexp(-(t2 -
tb1(i-2))**2/(0.5d0*
ta1(i-2))**2)
44 tl3 = dexp(-(t2 -
tb1(i-1))**2/(0.5d0*
ta1(i-1))**2)
45 tl4 = dexp(-(t2 -
tb1(i))**2/(0.5d0*
ta1(i))**2)
46 tl5 = dexp(-(t2 -
tb1(i+1))**2/(0.5d0*
ta1(i+1))**2)
47 tl6 = dexp(-(t2 -
tb1(i+2))**2/(0.5d0*
ta1(i+2))**2)
48 {
do ix^db = ixomin^db,ixomax^db\}
49 bp = floor((x(ix^d,1) - xprobmin1) / stt) + 1
51 if(x(ix^d,1) .lt. 0.d0)
then
52 lqgrid(ix^d) = tl1 * (
va1(bp + (i-4)*vlim)*(one-va0)+va0) &
53 + tl2 * (
va1(bp + (i-3)*vlim)*(one-va0)+va0) &
54 + tl3 * (
va1(bp + (i-2)*vlim)*(one-va0)+va0) &
55 + tl4 * (
va1(bp + (i-1)*vlim)*(one-va0)+va0) &
56 + tl5 * (
va1(bp + (i)*vlim)*(one-va0)+va0) &
57 + tl6 * (
va1(bp + (i+1)*vlim)*(one-va0)+va0)
65 if(qt .lt. trelax +
tb2(i))
then
66 tr1 = dexp(-(t2 -
tb2(i-3))**2/(0.5d0*
ta2(i-3))**2)
67 tr2 = dexp(-(t2 -
tb2(i-2))**2/(0.5d0*
ta2(i-2))**2)
68 tr3 = dexp(-(t2 -
tb2(i-1))**2/(0.5d0*
ta2(i-1))**2)
69 tr4 = dexp(-(t2 -
tb2(i))**2/(0.5d0*
ta2(i))**2)
70 tr5 = dexp(-(t2 -
tb2(i+1))**2/(0.5d0*
ta2(i+1))**2)
71 tr6 = dexp(-(t2 -
tb2(i+2))**2/(0.5d0*
ta2(i+2))**2)
72 {
do ix^db = ixomin^db,ixomax^db\}
73 bp = floor((x(ix^d,1) - xprobmin1) / stt) + 1
75 if(x(ix^d,1) .ge. 0.d0)
then
76 lqgrid(ix^d) = tr1 * (
va2(bp + (i-4)*vlim)*(one-va0)+va0) &
77 + tr2 * (
va2(bp + (i-3)*vlim)*(one-va0)+va0) &
78 + tr3 * (
va2(bp + (i-2)*vlim)*(one-va0)+va0) &
79 + tr4 * (
va2(bp + (i-1)*vlim)*(one-va0)+va0) &
80 + tr5 * (
va2(bp + (i)*vlim)*(one-va0)+va0) &
81 + tr6 * (
va2(bp + (i+1)*vlim)*(one-va0)+va0)
89 if(qt .lt. trelax)
then
93 else if(qt .lt. trelax +
tb1(1))
then
94 tr1 = dsin(t2 * dpi /
tb1(2))
97 else if(qt .lt. trelax +
tb1(2))
then
98 tr1 = dsin(t2 * dpi /
tb1(2))
99 tr3 = dsin((t2 -
tb1(1)) * dpi / (
tb1(3) -
tb1(1)))
103 if(t2 .lt.
tb1(i))
then
108 tr1 = dsin((t2 -
tb1(bt1-1)) * dpi / (
tb1(bt1+1) -
tb1(bt1-1)))
109 tr3 = dsin((t2 -
tb1(bt1)) * dpi / (
tb1(bt1+2) -
tb1(bt1)))
111 if(qt .lt. trelax)
then
115 else if(qt .lt. trelax +
tb2(1))
then
116 tr2 = dsin(t2 * dpi /
tb2(2))
119 else if(qt .lt. trelax +
tb2(2))
then
120 tr1 = dsin(t2 * dpi /
tb2(2))
121 tr3 = dsin((t2 -
tb2(1)) * dpi / (
tb2(3) -
tb2(1)))
125 if(t2 .lt.
tb2(i))
then
130 tr2 = dsin((t2 -
tb2(bt2-1)) * dpi / (
tb2(bt2+1) -
tb2(bt2-1)))
131 tr4 = dsin((t2 -
tb2(bt2)) * dpi / (
tb2(bt2+2) -
tb2(bt2)))
133 {
do ix^db = ixomin^db,ixomax^db\}
134 bp = floor((x(ix^d,1) - xprobmin1) / stt) + 1
135 if(x(ix^d,1) .lt. 0.d0)
then
136 lqgrid(ix^d) = tr1 * (
va1(bp + bt1*vlim)*(one-va0)+va0) + tr3 * (
va1(bp+(bt1+1)*vlim)*(one-va0)+va0)
138 lqgrid(ix^d) = tr2 * (
va2(bp + bt2*vlim)*(one-va0)+va0) + tr4 * (
va2(bp+(bt2+1)*vlim)*(one-va0)+va0)
147 subroutine generatetv(vtim,vlim,periods,variation,mnx,num,si)
148 integer,
intent(in) :: vtim,vlim,mnx,num
149 double precision,
intent(in) :: periods,variation,si
151 character(len=100) :: filename
161 if (
ndim/=2)
call mpistop(
"Randomized heating for 2D, 2.5D simulations only!")
163 write(filename,
"(a)")
"LeftT.dat"
164 inquire(file=filename,exist=alive)
166 open(unit=21,file=filename,status=
'old')
172 call randomt(
tb1,vtim,periods,variation,filename)
175 write(filename,
"(a)")
"RightT.dat"
176 inquire(file=filename,exist=alive)
178 open(unit=21,file=filename,status=
'old')
184 call randomt(
tb2,vtim,periods,variation,filename)
208 write(filename,
"(a)")
"LeftV.dat"
209 inquire(file=filename,exist=alive)
211 open(unit=22,file=filename,status=
'old')
217 call randomv(
va1,vtim,vlim,mnx,num,si,filename)
220 write(filename,
"(a)")
"RightV.dat"
221 inquire(file=filename,exist=alive)
223 open(unit=22,file=filename,status=
'old')
229 call randomv(
va2,vtim,vlim,mnx,num,si,filename)
240 subroutine randomt(tb,vtim,periods,variation,filename)
241 integer,
intent(in) :: vtim
242 double precision,
intent(in) :: periods,variation
243 double precision,
intent(inout) :: tb(:)
244 character(len=100),
intent(in) :: filename
246 double precision :: tt1,tt,mm
252 call random_number(mm)
253 tt1 = periods + variation * 2. * (mm - 0.5)
258 open(unit=21,file=filename,form=
'formatted',status=
'new')
259 write(21,
'(es12.4)') tb
264 subroutine randomv(va,vtim,vlim,mnx,num,si,filename)
265 integer,
intent(in) :: vtim,vlim,mnx,num
266 double precision,
intent(in) :: si
267 double precision,
intent(inout) :: va(:)
268 character(len=100),
intent(in) :: filename
269 double precision,
allocatable :: vn(:,:),vm(:)
270 double precision :: lambda,c,lambda_min,lambda_max,rms,ll
272 integer :: j, kk, lth = 1
274 allocate(vn(num, mnx))
277 lambda_min = lth * 1.d0 / (mnx * 1.d0)
278 lambda_max = lth * 2.d0
280 lambda = i * 1.d0 / num * (lambda_max - lambda_min) + lambda_min
282 call random_number(ll)
284 vn(i,j) = c*sin(2.*dpi*(j*1.d0/(mnx-1)*lth)/lambda+2.*dpi*2.*(ll-0.5))
291 vm(j) = vm(j) + vn(i, j)
295 vm = vm / dsqrt(rms/mnx)
297 va(vlim*(kk-1)+j) = vm(j)**2
301 open(unit=22,file=filename,form=
'formatted',status=
'new')
302 write(22,
'(es12.4)') va
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...
double precision unit_time
Physical scaling factor for time.
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 va2
double precision, dimension(:), allocatable tb1
double precision, dimension(:), allocatable ta2
double precision, dimension(:), allocatable ta1
double precision, dimension(:), allocatable va1
subroutine generatetv(vtim, vlim, periods, variation, mnx, num, si)
subroutine getlqgrid(lQgrid, ixIL, ixOL, qt, trelax, x, mode, vtim, vlim)
This module is for imposing randomized heating at x-axis for 2.5D simulation, for details please see ...
subroutine randomv(va, vtim, vlim, mnx, num, si, filename)
double precision, dimension(:), allocatable tb2
subroutine randomt(tb, vtim, periods, variation, filename)