23 double precision,
public ::
fld_mu = 0.6d0
66 integer,
allocatable,
public ::
i_opf(:)
69 double precision,
private,
protected :: fld_gamma
95 character(len=*),
intent(in) :: files(:)
104 do n = 1,
size(files)
105 open(
unitpar, file=trim(files(n)), status=
"old")
106 read(
unitpar, fld_list,
end=111)
119 subroutine fld_init(He_abundance, rhd_radiation_diffusion, rhd_gamma)
126 double precision,
intent(in) :: he_abundance, rhd_gamma
127 logical,
intent(in) :: rhd_radiation_diffusion
128 double precision :: sigma_thomson
131 character(len=1) :: ind_1
132 character(len=1) :: ind_2
133 character(len=2) :: cmp_f
134 character(len=5) :: cmp_e
143 write(ind_1,
'(I1)') idir
149 if (rhd_radiation_diffusion)
then
158 mg%operator_type = mg_vhelmholtz
166 fld_mu = (1.+4*he_abundance)/(2.+3.*he_abundance)
169 fld_gamma = rhd_gamma
174 sigma_thomson = 6.6524585
d-25
175 fld_kappa0 = sigma_thomson/const_mp * (1.+2.*he_abundance)/(1.+4.*he_abundance)
182 energy,qsourcesplit,active)
188 integer,
intent(in) :: ixi^
l, ixo^
l
189 double precision,
intent(in) :: qdt, x(ixi^s,1:
ndim)
190 double precision,
intent(in) :: wct(ixi^s,1:nw)
191 double precision,
intent(inout) :: w(ixi^s,1:nw)
192 logical,
intent(in) :: energy,qsourcesplit
193 logical,
intent(inout) :: active
195 double precision :: radiation_forcect(ixo^s,1:
ndim)
196 double precision :: kappact(ixo^s)
197 double precision :: rad_fluxct(ixo^s,1:
ndim)
199 double precision :: div_v(ixi^s,1:
ndim,1:
ndim)
200 double precision :: edd(ixo^s,1:
ndim,1:
ndim)
201 double precision :: nabla_vp(ixo^s)
202 double precision :: vel(ixi^s), grad_v(ixi^s), grad0_v(ixo^s)
203 double precision :: grad_e(ixo^s)
205 integer :: idir, jdir
207 double precision :: fld_r(ixo^s), lambda(ixo^s)
225 w(ixo^s,iw_mom(idir)) = w(ixo^s,iw_mom(idir)) &
226 + qdt * wct(ixo^s,iw_rho)*radiation_forcect(ixo^s,idir)
232 w(ixo^s,iw_e) = w(ixo^s,iw_e) &
233 + qdt * wct(ixo^s,iw_mom(idir))*radiation_forcect(ixo^s,idir)
244 vel(ixi^s) = wct(ixi^s,iw_mom(jdir))/wct(ixi^s,iw_rho)
247 div_v(ixo^s,idir,jdir) = grad_v(ixo^s)
259 nabla_vp(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1) &
264 nabla_vp(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1) &
265 + div_v(ixo^s,1,2)*edd(ixo^s,1,2) &
266 + div_v(ixo^s,2,1)*edd(ixo^s,2,1) &
267 + div_v(ixo^s,2,2)*edd(ixo^s,2,2)
271 nabla_vp(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1) &
272 + div_v(ixo^s,1,2)*edd(ixo^s,1,2) &
273 + div_v(ixo^s,1,3)*edd(ixo^s,1,3) &
274 + div_v(ixo^s,2,1)*edd(ixo^s,2,1) &
275 + div_v(ixo^s,2,2)*edd(ixo^s,2,2) &
276 + div_v(ixo^s,2,3)*edd(ixo^s,2,3) &
277 + div_v(ixo^s,3,1)*edd(ixo^s,3,1) &
278 + div_v(ixo^s,3,2)*edd(ixo^s,3,2) &
279 + div_v(ixo^s,3,3)*edd(ixo^s,3,3)
282 w(ixo^s,iw_r_e) = w(ixo^s,iw_r_e) &
283 - qdt * nabla_vp(ixo^s)*wct(ixo^s,iw_r_e)
294 integer,
intent(in) :: ixi^
l, ixo^
l
295 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim), w(ixi^s,1:nw)
296 double precision,
intent(inout) :: dtnew
298 double precision :: radiation_force(ixo^s,1:
ndim)
299 double precision :: rad_flux(ixo^s,1:
ndim)
300 double precision :: kappa(ixo^s)
301 double precision :: dxinv(1:
ndim), max_grad
304 ^
d&dxinv(^
d)=one/
dx^
d;
310 radiation_force(ixo^s,idim) = w(ixo^s,iw_rho)*kappa(ixo^s)*rad_flux(ixo^s,idim)/(const_c/
unit_velocity)
311 max_grad = maxval(abs(radiation_force(ixo^s,idim)))
312 max_grad = max(max_grad, epsilon(1.0d0))
313 dtnew = min(dtnew,
courantpar / sqrt(max_grad * dxinv(idim)))
321 energy,qsourcesplit,active)
328 integer,
intent(in) :: ixi^
l, ixo^
l
329 double precision,
intent(in) :: qdt, x(ixi^s,1:
ndim)
330 double precision,
intent(in) :: wct(ixi^s,1:nw)
331 double precision,
intent(inout) :: w(ixi^s,1:nw)
332 logical,
intent(in) :: energy,qsourcesplit
333 logical,
intent(inout) :: active
352 integer,
intent(in) :: ixi^
l, ixo^
l
353 double precision,
intent(in) :: w(ixi^s, 1:nw)
354 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
355 double precision,
intent(out) :: fld_kappa(ixo^s)
356 double precision :: temp(ixi^s), pth(ixi^s), a2(ixo^s)
357 double precision :: rho0,temp0,n,sigma_b
358 double precision :: akram, bkram
359 double precision :: vth(ixo^s), gradv(ixi^s), eta(ixo^s), t(ixo^s)
361 integer :: i,j,ix^
d, idir
379 call phys_get_pthermal(w,x,ixi^
l,ixo^
l,temp)
380 temp(ixo^s)=temp(ixo^s)/w(ixo^s,iw_rho)
387 call phys_get_pthermal(w,x,ixi^
l,ixo^
l,pth)
390 akram = 13.1351597305
391 bkram = -4.5182188206
394 * (1.d0+10.d0**akram*w(ixo^s,iw_rho)*
unit_density*(a2(ixo^s)/1.d12)**bkram)
396 {
do ix^
d=ixomin^
d,ixomax^
d\ }
408 {
do ix^
d=ixomin^
d,ixomax^
d\ }
417 call mpistop(
"special opacity not defined")
422 call mpistop(
"Doesn't know opacity law")
435 integer,
intent(in) :: ixi^
l, ixo^
l
436 double precision,
intent(in) :: w(ixi^s, 1:nw)
437 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
438 double precision,
intent(out) :: fld_r(ixo^s), fld_lambda(ixo^s)
439 double precision :: kappa(ixo^s)
440 double precision :: normgrad2(ixo^s)
441 double precision :: grad_r_e(ixi^s), grad_r_eo(ixo^s)
442 double precision :: rad_e(ixi^s)
443 integer :: idir, ix^
d
445 double precision :: tmp_l(ixi^s), filtered_l(ixi^s)
446 integer :: filter, idim
450 fld_lambda(ixo^s) = 1.d0/3.d0
456 normgrad2(ixo^s) = 0.d0
458 rad_e(ixi^s) = w(ixi^s, iw_r_e)
461 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_eo(ixo^s)**2
469 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
472 fld_lambda(ixo^s) = one/fld_r(ixo^s)
477 normgrad2(ixo^s) = 0.d0
479 rad_e(ixi^s) = w(ixi^s, iw_r_e)
482 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_eo(ixo^s)**2
490 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
494 fld_lambda(ixo^s) = (2.d0+fld_r(ixo^s))/(6.d0+3*fld_r(ixo^s)+fld_r(ixo^s)**2.d0)
497 call mpistop(
"Pomraning2 is not quite working, use Pomraning or Minerbo")
500 normgrad2(ixo^s) = 0.d0
502 rad_e(ixi^s) = w(ixi^s, iw_r_e)
505 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_eo(ixo^s)**2
513 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
517 fld_lambda(ixo^s) = one/fld_r(ixo^s)*(one/dtanh(fld_r(ixo^s)) - one/fld_r(ixo^s))
521 where(dtanh(fld_r(ixo^s)) .ne. dtanh(fld_r(ixo^s)))
522 fld_lambda(ixo^s) = 1.d0/3.d0
528 normgrad2(ixo^s) = 0.d0
530 rad_e(ixi^s) = w(ixi^s, iw_r_e)
533 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_eo(ixo^s)**2
541 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
545 {
do ix^
d=ixomin^
d,ixomax^
d\ }
546 if (fld_r(ix^
d) .lt. 3.d0/2.d0)
then
547 fld_lambda(ix^
d) = 2.d0/(3.d0 + dsqrt(9.d0 + 12.d0*fld_r(ix^
d)**2.d0))
549 fld_lambda(ix^
d) = 1.d0/(1.d0 + fld_r(ix^
d) + dsqrt(1.d0 + 2.d0*fld_r(ix^
d)))
555 call mpistop(
"special fluxlimiter not defined")
559 call mpistop(
'Fluxlimiter unknown')
564 if (
size_l_filter .lt. 1)
call mpistop(
"D filter of size < 1 makes no sense")
567 tmp_l(ixo^s) = fld_lambda(ixo^s)
568 filtered_l(ixo^s) = zero
573 filtered_l(ix^
d) = filtered_l(ix^
d) &
574 + tmp_l(ix^
d+filter*
kr(idim,^
d)) &
575 + tmp_l(ix^
d-filter*
kr(idim,^
d))
584 fld_lambda(ixo^s) = tmp_l(ixo^s)
595 integer,
intent(in) :: ixi^
l, ixo^
l
596 double precision,
intent(in) :: w(ixi^s, 1:nw)
597 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
598 double precision,
intent(out) :: rad_flux(ixo^s, 1:
ndim)
600 double precision :: grad_r_e(ixi^s), grad_r_eo(ixo^s)
601 double precision :: rad_e(ixi^s)
602 double precision :: kappa(ixo^s), lambda(ixo^s), fld_r(ixo^s)
603 integer :: ix^
d, idir
605 rad_e(ixi^s) = w(ixi^s, iw_r_e)
614 rad_flux(ixo^s, idir) = -(const_c/
unit_velocity)*lambda(ixo^s)/(kappa(ixo^s)*w(ixo^s,iw_rho))*grad_r_eo(ixo^s)
628 integer,
intent(in) :: ixI^L, ixO^L
629 double precision,
intent(in) :: w(ixI^S, 1:nw)
630 double precision,
intent(in) :: x(ixI^S, 1:ndim)
631 double precision,
intent(out) :: eddington_tensor(ixO^S,1:ndim,1:ndim)
632 double precision :: tnsr2(ixO^S,1:ndim,1:ndim)
633 double precision :: normgrad2(ixO^S), f(ixO^S)
634 double precision :: grad_r_e(ixI^S, 1:ndim), rad_e(ixI^S)
635 double precision :: grad_r_eO(ixO^S, 1:ndim)
636 double precision :: lambda(ixO^S), fld_R(ixO^S)
637 integer :: i,j, idir,jdir
641 normgrad2(ixo^s) = zero
643 rad_e(ixi^s) = w(ixi^s, iw_r_e)
646 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_eo(ixo^s,idir)**2
656 f(ixo^s) = lambda(ixo^s) + lambda(ixo^s)**two * fld_r(ixo^s)**two
657 f(ixo^s) = one/two*(one-f(ixo^s)) + one/two*(3.d0*f(ixo^s) - one)
660 eddington_tensor(ixo^s,1,1) = f(ixo^s)
665 eddington_tensor(ixo^s,idir,idir) = half*(one-f(ixo^s))
670 if (idir .ne. jdir) eddington_tensor(ixo^s,idir,jdir) = zero
671 tnsr2(ixo^s,idir,jdir) = half*(3.d0*f(ixo^s) - 1)&
672 *grad_r_eo(ixo^s,idir)*grad_r_eo(ixo^s,jdir)/normgrad2(ixo^s)
686 where ((tnsr2(ixo^s,idir,jdir) .eq. tnsr2(ixo^s,idir,jdir)) &
687 .and. (normgrad2(ixo^s) .gt. smalldouble))
688 eddington_tensor(ixo^s,idir,jdir) = eddington_tensor(ixo^s,idir,jdir) + tnsr2(ixo^s,idir,jdir)
701 integer,
intent(in) :: ixi^
l, ixo^
l
702 double precision,
intent(in) :: w(ixi^s, 1:nw)
703 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
704 double precision :: eddington_tensor(ixo^s,1:
ndim,1:
ndim)
705 double precision,
intent(out):: rad_pressure(ixo^s,1:
ndim,1:
ndim)
713 rad_pressure(ixo^s,i,j) = eddington_tensor(ixo^s,i,j)* w(ixo^s,iw_r_e)
731 type(state),
target :: psa(max_blocks)
732 type(state),
target :: psb(max_blocks)
733 double precision,
intent(in) :: qdt
734 double precision,
intent(in) :: qtC
735 double precision,
intent(in) :: dtfactor
738 double precision :: res, max_residual, lambda
739 integer,
parameter :: max_its = 100
741 integer :: iw_to,iw_from
742 integer :: iigrid, igrid, id
743 integer :: nc, lvl, i
745 real(dp) :: fac, facD
755 if (qdt <
dtmin)
then
757 print *,
'skipping implicit solve: dt too small!'
765 mg%operator_type = mg_vhelmholtz
766 mg%smoother_type = mg_smoother_gs
767 call mg_set_methods(
mg)
769 if (.not.
mg%is_allocated)
call mpistop(
"multigrid tree not allocated yet")
772 lambda = 1.d0/(dtfactor *
dt_diff)
773 call vhelmholtz_set_lambda(lambda)
877 call mg_restrict(
mg, mg_iveps)
878 call mg_fill_ghost_cells(
mg, mg_iveps)
884 call mg_fas_fmg(
mg, .true., max_res=res)
887 if (res < max_residual)
exit
888 call mg_fas_vcycle(
mg, max_res=res)
894 if (res .le. 0.d0)
then
896 if (
mg%my_rank == 0) &
897 write(*,*)
it,
' resiudal zero ', res
900 if (
mg%my_rank == 0)
then
902 error stop
"Diffusion residual to zero"
909 if (n == max_its + 1)
then
911 if (
mg%my_rank == 0) &
912 write(*,*)
it,
' resiudal high ', res
915 if (
mg%my_rank == 0)
then
916 print *,
"Did you specify boundary conditions correctly?"
917 print *,
"Or is the variation in diffusion too large?"
919 print *,
mg%bc(1, mg_iphi)%bc_value,
mg%bc(2, mg_iphi)%bc_value
921 error stop
"diffusion_solve: no convergence"
970 type(state),
target :: psa(max_blocks)
971 double precision,
intent(in) :: qtC
973 integer :: iigrid, igrid, level
980 do iigrid=1,igridstail; igrid=igrids(iigrid);
995 integer,
intent(in) :: ixI^L, ixO^L
996 double precision,
intent(inout) :: w(ixI^S, 1:nw)
998 double precision :: gradE(ixO^S),divF(ixO^S)
999 double precision :: divF_h(ixO^S),divF_j(ixO^S)
1000 double precision :: diff_term(ixO^S)
1002 integer :: idir, jxO^L, hxO^L
1008 hxo^l=ixo^l-
kr(idir,^
d);
1009 jxo^l=ixo^l+
kr(idir,^
d);
1013 divf(ixo^s) = divf(ixo^s) + 2.d0*(divf_h(ixo^s) + divf_j(ixo^s))/
dxlevel(idir)**2
1016 w(ixo^s,iw_r_e) = divf(ixo^s)
1027 integer,
intent(in) :: ixI^L, ixO^L
1028 double precision,
intent(inout) :: w(ixI^S, 1:nw)
1029 double precision,
intent(in) :: wCT(ixI^S, 1:nw)
1030 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1032 double precision :: kappa(ixO^S), lambda(ixO^S), fld_R(ixO^S)
1034 double precision :: max_D(ixI^S), grad_r_e(ixI^S), rad_e(ixI^S)
1035 integer :: idir,i,j, ix^D
1060 type(state),
target :: psa(max_blocks)
1063 integer :: iigrid, igrid, level
1068 do iigrid=1,igridstail; igrid=igrids(iigrid);
1082 integer,
intent(in) :: ixI^L, ixO^L
1083 double precision,
intent(inout) :: w(ixI^S, 1:nw)
1085 double precision :: tmp_D(ixI^S), filtered_D(ixI^S)
1086 integer :: ix^D, filter, idim
1088 if (
size_d_filter .lt. 1)
call mpistop(
"D filter of size < 1 makes no sense")
1092 filtered_d(ixo^s) = zero
1097 filtered_d(ix^d) = filtered_d(ix^d) &
1098 + tmp_d(ix^d+filter*
kr(idim,^d)) &
1099 + tmp_d(ix^d-filter*
kr(idim,^d))
1127 integer,
intent(in) :: ixI^L, ixO^L
1128 double precision,
intent(in) :: x(ixI^S,1:ndim)
1129 double precision,
intent(in) :: wCT(ixI^S,1:nw)
1130 double precision,
intent(inout) :: w(ixI^S,1:nw)
1132 double precision :: a1(ixO^S), a2(ixO^S)
1133 double precision :: c0(ixO^S), c1(ixO^S)
1134 double precision :: e_gas(ixO^S), E_rad(ixO^S)
1135 double precision :: kappa(ixO^S)
1136 double precision :: sigma_b
1138 integer :: i,j,idir,ix^D
1142 e_gas(ixo^s) = wct(ixo^s,iw_e) - half*sum(wct(ixo^s, iw_mom(:))**2, dim=ndim+1)/wct(ixo^s, iw_rho)
1147 {
do ix^d=ixomin^d,ixomax^d\ }
1151 e_rad(ixo^s) = wct(ixo^s,iw_r_e)
1162 a1(ixo^s) = const_rad_a*(
fld_mu*const_mp/const_kb*(fld_gamma-1))**4/(wct(ixo^s,iw_rho)*
unit_density)**4 &
1164 a2(ixo^s) = e_gas(ixo^s) + e_rad(ixo^s)
1166 c0(ixo^s) = a2(ixo^s)/a1(ixo^s)
1167 c1(ixo^s) = 1.d0/a1(ixo^s)
1170 a1(ixo^s) = 4*kappa(ixo^s)*sigma_b*(fld_gamma-one)**4/wct(ixo^s,iw_rho)**3*
dt
1173 c0(ixo^s) = ((one + a2(ixo^s))*e_gas(ixo^s) + a2(ixo^s)*e_rad(ixo^s))/a1(ixo^s)
1174 c1(ixo^s) = (one + a2(ixo^s))/a1(ixo^s)
1178 {
do ix^d=ixomin^d,ixomax^d\ }
1183 call newton_method(e_gas(ix^d), e_rad(ix^d), c0(ix^d), c1(ix^d))
1185 call halley_method(e_gas(ix^d), e_rad(ix^d), c0(ix^d), c1(ix^d))
1187 call halley_method(e_gas(ix^d), e_rad(ix^d), c0(ix^d), c1(ix^d))
1189 call mpistop(
'root-method not known')
1194 e_rad(ixo^s) = a2(ixo^s) - e_gas(ixo^s)
1197 e_rad(ixo^s) = (a1(ixo^s)*e_gas(ixo^s)**4.d0 + e_rad(ixo^s))/(one + a2(ixo^s))
1221 w(ixo^s,iw_e) = e_gas(ixo^s)
1226 w(ixo^s,iw_e) = w(ixo^s,iw_e) + half*sum(wct(ixo^s, iw_mom(:))**2, dim=ndim+1)/wct(ixo^s, iw_rho)
1233 w(ixo^s,iw_r_e) = e_rad(ixo^s)
1242 double precision,
intent(in) :: c0, c1
1243 double precision,
intent(in) :: E_rad
1244 double precision,
intent(inout) :: e_gas
1246 double precision :: bisect_a, bisect_b, bisect_c
1247 integer :: n, max_its
1253 bisect_b = min(abs(c0/c1),abs(c0)**(1.d0/4.d0))+smalldouble
1257 do while (abs(bisect_b-bisect_a) .ge.
fld_bisect_tol*min(e_gas,e_rad))
1258 bisect_c = (bisect_a + bisect_b)/two
1261 if (n .gt. max_its)
then
1263 call mpistop(
'No convergece in bisection scheme')
1288 call mpistop(
"Problem with fld bisection method")
1296 print*,
"IGNORING GAS-RAD ENERGY EXCHANGE ", c0, c1
1311 2435 e_gas = (bisect_a + bisect_b)/two
1318 double precision,
intent(in) :: c0, c1
1319 double precision,
intent(in) :: E_rad
1320 double precision,
intent(inout) :: e_gas
1322 double precision :: xval, yval, der, deltax
1337 xval = xval + deltax
1339 if (ii .gt. 1d3)
then
1340 print*,
'skip to bisection algorithm'
1353 double precision,
intent(in) :: c0, c1
1354 double precision,
intent(in) :: E_rad
1355 double precision,
intent(inout) :: e_gas
1357 double precision :: xval, yval, der, dder, deltax
1373 deltax = -two*yval*der/(two*der**2 - yval*dder)
1374 xval = xval + deltax
1376 if (ii .gt. 1d3)
then
1390 double precision,
intent(in) :: e_gas
1391 double precision,
intent(in) :: c0, c1
1392 double precision :: val
1394 val = e_gas**4.d0 + c1*e_gas - c0
1401 double precision,
intent(in) :: e_gas
1402 double precision,
intent(in) :: c0, c1
1403 double precision :: der
1405 der = 4.d0*e_gas**3.d0 + c1
1412 double precision,
intent(in) :: e_gas
1413 double precision,
intent(in) :: c0, c1
1414 double precision :: dder
1416 dder = 4.d0*3.d0*e_gas**2.d0
1424 integer,
intent(in) :: ixI^L, ixO^L, idir
1425 integer,
intent(in) :: n
1427 double precision,
intent(in) :: q(ixI^S), x(ixI^S,1:ndim)
1428 double precision,
intent(out) :: gradq(ixO^S)
1429 integer :: jxO^L, hxO^L
1443 call mpistop(
"gradientO stencil too wide")
1444 elseif (n .eq. 1)
then
1445 hxo^l=ixo^l-
kr(idir,^
d);
1446 jxo^l=ixo^l+
kr(idir,^
d);
1447 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
1448 elseif (n .eq. 2)
then
1451 hxo^l=ixo^l-
kr(idir,^
d);
1452 jxo^l=ixo^l+
kr(idir,^
d);
1453 gradq(ixo^s) = gradq(ixo^s) + 2.d0/3.d0*(q(jxo^s)-q(hxo^s))
1455 hxo^l=ixo^l-2*
kr(idir,^
d);
1456 jxo^l=ixo^l+2*
kr(idir,^
d);
1457 gradq(ixo^s) = gradq(ixo^s) - 1.d0/12.d0*(q(jxo^s)-q(hxo^s))
1459 gradq(ixo^s) = gradq(ixo^s)/
dxlevel(idir)
1460 elseif (n .eq. 3)
then
1463 hxo^l=ixo^l-
kr(idir,^
d);
1464 jxo^l=ixo^l+
kr(idir,^
d);
1465 gradq(ixo^s) = gradq(ixo^s) + 3.d0/4.d0*(q(jxo^s)-q(hxo^s))
1467 hxo^l=ixo^l-2*
kr(idir,^
d);
1468 jxo^l=ixo^l+2*
kr(idir,^
d);
1469 gradq(ixo^s) = gradq(ixo^s) - 3.d0/20.d0*(q(jxo^s)-q(hxo^s))
1471 hxo^l=ixo^l-3*
kr(idir,^
d);
1472 jxo^l=ixo^l+3*
kr(idir,^
d);
1473 gradq(ixo^s) = gradq(ixo^s) + 1.d0/60.d0*(q(jxo^s)-q(hxo^s))
1475 gradq(ixo^s) = gradq(ixo^s)/
dxlevel(idir)
1477 call mpistop(
"gradientO stencil unknown")
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter const_c
Nicolas Moens Module for including flux limited diffusion (FLD)-approximation in Radiation-hydrodynam...
subroutine energy_interaction(w, wCT, x, ixIL, ixOL)
This subroutine calculates the radiation heating, radiation cooling and photon tiring using an implic...
subroutine gradiento(q, x, ixIL, ixOL, idir, gradq, n)
Calculate gradient of a scalar q within ixL in direction idir difference with gradient is gradq(ixO^S...
double precision, public fld_mu
mean particle mass
double precision, public fld_bisect_tol
Tolerance for bisection method for Energy sourceterms This is a percentage of the minimum of gas- and...
double precision, public fld_diff_tol
Tolerance for adi method for radiative Energy diffusion.
subroutine diffuse_e_rad_mg(dtfactor, qdt, qtC, psa, psb)
Calling all subroutines to perform the multigrid method Communicates rad_e and diff_coeff to multigri...
subroutine fld_smooth_diffcoef(w, ixIL, ixOL)
Use running average on Diffusion coefficient.
subroutine, public fld_radforce_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine fld_params_read(files)
public methods these are called in mod_rhd_phys
logical diffcrash_resume
Resume run when multigrid returns error.
subroutine update_diffcoeff(psa)
integer, dimension(:), allocatable, public i_opf
Index for Flux weighted opacities.
character(len=8) fld_opacity_law
Use constant Opacity?
subroutine fld_get_eddington(w, x, ixIL, ixOL, eddington_tensor)
Calculate Eddington-tensor Stores Eddington-tensor in w-array.
double precision function ddpolynomial_bisection(e_gas, c0, c1)
Evaluate second derivative of polynomial at argument e_gas.
subroutine, public fld_get_radflux(w, x, ixIL, ixOL, rad_flux)
Calculate Radiation Flux stores radiation flux in w-array.
character(len=16) fld_fluxlimiter
Diffusion limit lambda = 0.33.
subroutine halley_method(e_gas, E_rad, c0, c1)
Find the root of the 4th degree polynomial using the Halley method.
subroutine put_diffterm_onegrid(ixIL, ixOL, w)
inplace update of psa==>F_im(psa)
character(len=50) fld_opal_table
double precision, public fld_kappa0
Opacity per unit of unit_density.
character(len=8) fld_diff_scheme
Which method to solve diffusion part.
logical lineforce_opacities
Use or dont use lineforce opacities.
subroutine, public get_fld_rad_force(qdt, ixIL, ixOL, wCT, w, x, energy, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO This subroutine handles th...
double precision function polynomial_bisection(e_gas, c0, c1)
Evaluate polynomial at argument e_gas.
logical fld_eint_split
source split for energy interact and radforce:
character(len=8) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
subroutine, public fld_get_radpress(w, x, ixIL, ixOL, rad_pressure)
Calculate Radiation Pressure Returns Radiation Pressure as tensor.
logical fld_radforce_split
subroutine, public fld_init(He_abundance, rhd_radiation_diffusion, rhd_gamma)
Initialising FLD-module: Read opacities Initialise Multigrid adimensionalise kappa Add extra variable...
double precision function dpolynomial_bisection(e_gas, c0, c1)
Evaluate first derivative of polynomial at argument e_gas.
subroutine, public fld_get_opacity(w, x, ixIL, ixOL, fld_kappa)
Sets the opacity in the w-array by calling mod_opal_opacity.
subroutine, public get_fld_energy_interact(qdt, ixIL, ixOL, wCT, w, x, energy, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO This subroutine handles th...
integer i_diff_mg
diffusion coefficient for multigrid method
subroutine evaluate_e_rad_mg(qtC, psa)
inplace update of psa==>F_im(psa)
subroutine fld_get_diffcoef_central(w, wCT, x, ixIL, ixOL)
Calculates cell-centered diffusion coefficient to be used in multigrid.
subroutine bisection_method(e_gas, E_rad, c0, c1)
Find the root of the 4th degree polynomial using the bisection method.
logical flux_lim_filter
Take a running average over the fluxlimiter.
logical diff_coef_filter
Take running average for Diffusion coefficient.
subroutine, public fld_get_fluxlimiter(w, x, ixIL, ixOL, fld_lambda, fld_R)
Calculate fld flux limiter This subroutine calculates flux limiter lambda using the prescription stor...
double precision dt_diff
running timestep for diffusion solver, initialised as zero
double precision, public diff_crit
Number for splitting the diffusion module.
subroutine newton_method(e_gas, E_rad, c0, c1)
Find the root of the 4th degree polynomial using the Newton-Ralphson method.
Module with basic grid data structures.
Module with geometry-related routines (e.g., divergence, curl)
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
update ghost cells of all blocks including physical boundaries
subroutine getbc(time, qdt, psb, nwstart, nwbc, req_diag)
do update ghost cells of all blocks including physical boundaries
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision small_pressure
double precision unit_density
Physical scaling factor for density.
double precision unit_opacity
Physical scaling factor for Opacity.
integer, parameter unitpar
file handle for IO
double precision global_time
The global simulation time.
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer it
Number of time steps taken.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter ndim
Number of spatial dimensions for grid variables.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
integer, dimension(:), allocatable, parameter d
double precision courantpar
The Courant (CFL) number used for the simulation.
double precision unit_velocity
Physical scaling factor for velocity.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision unit_temperature
Physical scaling factor for temperature.
double precision, dimension(:,:), allocatable dx
integer nghostcells
Number of ghost cells surrounding a grid.
logical use_multigrid
Use multigrid (only available in 2D and 3D)
double precision dtmin
Stop the simulation when the time step becomes smaller than this value.
double precision, dimension(ndim) dxlevel
Module to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
type(mg_t) mg
Data structure containing the multigrid tree.
subroutine mg_copy_to_tree(iw_from, iw_to, restrict, restrict_gc, factor, state_from)
Copy a variable to the multigrid tree, including a layer of ghost cells.
subroutine mg_copy_from_tree_gc(iw_from, iw_to, state_to)
Copy from multigrid tree with one layer of ghost cells. Corner ghost cells are not used/set.
This module reads in Rosseland-mean opacities from OPAL tables. Table opacity values are given in bas...
subroutine, public init_opal_table(tablename)
This subroutine is called when the FLD radiation module is initialised. The OPAL tables for different...
subroutine, public set_opal_opacity(rho, temp, kappa)
This subroutine calculates the opacity for a given temperature-density structure. Opacities are read ...
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_get_tgas), pointer phys_get_tgas
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
logical phys_req_diagonal
Whether the physics routines require diagonal ghost cells, for example for computing a curl.
procedure(sub_get_pthermal), pointer phys_get_pthermal
procedure(sub_implicit_update), pointer phys_implicit_update
Module with all the methods that users can customize in AMRVAC.
procedure(special_opacity), pointer usr_special_opacity
procedure(special_diffcoef), pointer usr_special_diffcoef
procedure(special_fluxlimiter), pointer usr_special_fluxlimiter
procedure(special_opacity_qdot), pointer usr_special_opacity_qdot
integer function var_set_extravar(name_cons, name_prim, ix)
Set extra variable in w, which is not advected and has no boundary conditions. This has to be done af...
The data structure that contains information about a tree node/grid block.