40 double precision,
private,
protected :: fld_gamma
60 character(len=*),
intent(in) :: files(:)
68 open(
unitpar, file=trim(files(n)), status=
"old")
69 read(
unitpar, fld_list,
end=111)
84 double precision,
intent(in) :: he_abundance, r_gamma
85 double precision :: sigma_thomson
99 call mpistop(
"nth_for_diff_mg must be 1 or 2")
109 mg%operator_type = mg_vhelmholtz
110 mg%smoother_type = mg_smoother_gs
118 sigma_thomson = 6.6524585
d-25
119 fld_kappa0 = sigma_thomson/const_mp * (1.+2.*he_abundance)/(1.+4.*he_abundance)
136 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
137 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
140 mg%bc(ib, mg_iphi)%bc_type = mg_bc_dirichlet
141 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
144 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
145 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
150 call mpistop(
"special BC for MG not defined")
156 call mpistop(
"special BC for MG not defined")
161 call mpistop(
"divE_multigrid warning: unknown b.c. ")
165 mg%bc(ib, mg_iveps)%bc_type = mg_bc_neumann
166 mg%bc(ib, mg_iveps)%bc_value = 0.0_dp
181 integer,
intent(in) :: ixi^
l, ixo^
l
182 double precision,
intent(in) :: qdt, x(ixi^s,1:
ndim)
183 double precision,
intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw)
184 double precision,
intent(inout) :: w(ixi^s,1:nw)
185 logical,
intent(in) :: qsourcesplit
186 logical,
intent(inout) :: active
188 integer :: idir,jdir,nth_for_fld,ix^
d
189 double precision :: sigma_b,cc
190 double precision,
dimension(ixI^S) :: a1,a2,a3,c0,c1,kappa
191 double precision,
dimension(ixI^S) :: e_gas,e_rad,tmp
192 double precision,
dimension(ixI^S,1:ndim,1:ndim) :: div_v,edd
201 call gradient(wctprim(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,tmp,nth_for_fld)
204 tmp(ixo^s) = -a1(ixo^s)*tmp(ixo^s)
206 w(ixo^s,iw_mom(idir)) = w(ixo^s,iw_mom(idir))+ qdt*tmp(ixo^s)
208 w(ixo^s,iw_e) = w(ixo^s,iw_e) + qdt*wctprim(ixo^s,iw_mom(idir))*tmp(ixo^s)
215 call gradient(wctprim(ixi^s,iw_mom(jdir)),ixi^
l,ixo^
l,idir,tmp)
216 div_v(ixo^s,idir,jdir) = tmp(ixo^s)
221 a3(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1)
224 a3(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1) &
225 + div_v(ixo^s,1,2)*edd(ixo^s,1,2) &
226 + div_v(ixo^s,2,1)*edd(ixo^s,2,1) &
227 + div_v(ixo^s,2,2)*edd(ixo^s,2,2)
230 a3(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1) &
231 + div_v(ixo^s,1,2)*edd(ixo^s,1,2) &
232 + div_v(ixo^s,1,3)*edd(ixo^s,1,3) &
233 + div_v(ixo^s,2,1)*edd(ixo^s,2,1) &
234 + div_v(ixo^s,2,2)*edd(ixo^s,2,2) &
235 + div_v(ixo^s,2,3)*edd(ixo^s,2,3) &
236 + div_v(ixo^s,3,1)*edd(ixo^s,3,1) &
237 + div_v(ixo^s,3,2)*edd(ixo^s,3,2) &
238 + div_v(ixo^s,3,3)*edd(ixo^s,3,3)
245 e_gas(ixo^s) = w(ixo^s,iw_e)-
half*sum(w(ixo^s,iw_mom(:))**2,dim=
ndim+1)/w(ixo^s,iw_rho)
246 if(
allocated(iw_mag))
then
247 e_gas(ixo^s) = e_gas(ixo^s)-
half*sum(w(ixo^s,iw_mag(:))**2,dim=
ndim+1)
249 e_rad(ixo^s) = w(ixo^s,iw_r_e)
252 {
do ix^
d = ixomin^
d,ixomax^
d\ }
260 a1(ixo^s) = 4.d0*sigma_b*qdt*kappa(ixo^s)*(fld_gamma-
one)**4/(wct(ixo^s,iw_rho)**3*tmp(ixo^s)**4)
261 a2(ixo^s) = cc*kappa(ixo^s)*wct(ixo^s,iw_rho)*qdt
262 a3(ixo^s) = a3(ixo^s)*qdt
264 c0(ixo^s) = ((
one+a2(ixo^s)+a3(ixo^s))*e_gas(ixo^s)+a2(ixo^s)*e_rad(ixo^s))/a1(ixo^s)/(
one+a3(ixo^s))
265 c1(ixo^s) = (
one+a2(ixo^s)+a3(ixo^s))/a1(ixo^s)/(
one+a3(ixo^s))
268 {
do ix^
d = ixomin^
d,ixomax^
d\}
277 call mpistop(
'root-method not known')
281 e_rad(ixo^s) = (a1(ixo^s)*e_gas(ixo^s)**4.d0+e_rad(ixo^s))/(
one+a2(ixo^s)+a3(ixo^s))
284 {
do ix^
d = ixomin^
d,ixomax^
d\ }
291 w(ixo^s,iw_e) = e_gas(ixo^s)
292 w(ixo^s,iw_e) = w(ixo^s,iw_e)+
half*sum(w(ixo^s,iw_mom(:))**2,dim=
ndim+1)/w(ixo^s,iw_rho)
293 if(
allocated(iw_mag))
then
294 w(ixo^s,iw_e) = w(ixo^s,iw_e)+
half*sum(w(ixo^s,iw_mag(:))**2,dim=
ndim+1)
297 w(ixo^s,iw_r_e) = e_rad(ixo^s)
309 integer,
intent(in) :: ixi^
l, ixo^
l
310 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim), w(ixi^s,1:nw)
311 double precision,
intent(inout) :: dtnew
313 integer :: idim,idims,nth_for_fld
314 double precision :: dxinv(1:
ndim), max_grav
315 double precision :: lambda(ixi^s),fld_r(ixi^s)
316 double precision :: tmp(ixi^s)
317 double precision :: cmax(ixi^s),cmaxtot(ixi^s),courantmaxtots
319 if(
fld_debug)print *,
'DT limit on entry to radforce_get_dt=',dtnew
323 ^
d&dxinv(^
d)=one/
dx^
d;
325 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idim,tmp,nth_for_fld)
326 max_grav = maxval(dabs(-lambda(ixo^s)*tmp(ixo^s)/w(ixo^s,iw_rho)))
327 max_grav = max(max_grav, epsilon(1.0d0))
328 dtnew = min(dtnew, 1.0d0 / dsqrt(max_grav * dxinv(idim)))
332 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idim,tmp,nth_for_fld)
333 max_grav = maxval(dabs(-lambda(ixo^s)*tmp(ixo^s)/w(ixo^s,iw_rho))/
block%ds(ixo^s,idim))
334 max_grav = max(max_grav, epsilon(1.0d0))
335 dtnew = min(dtnew, 1.0d0 / dsqrt(max_grav))
339 if(
fld_debug)print *,
'DT limit after RADFORCE eff grav=',dtnew
344 ^
d&dxinv(^
d)=one/
dx^
d;
346 cmax(ixo^s)=dabs(w(ixo^s,iw_mom(idims)))+dsqrt(tmp(ixo^s))
348 cmaxtot(ixo^s)=cmax(ixo^s)*dxinv(idims)
350 cmaxtot(ixo^s)=cmaxtot(ixo^s)+cmax(ixo^s)*dxinv(idims)
355 cmax(ixo^s)=dabs(w(ixo^s,iw_mom(idims)))+dsqrt(tmp(ixo^s))
357 cmaxtot(ixo^s)=cmax(ixo^s)/
block%ds(ixo^s,idims)
359 cmaxtot(ixo^s)=cmaxtot(ixo^s)+cmax(ixo^s)/
block%ds(ixo^s,idims)
364 courantmaxtots=maxval(cmaxtot(ixo^s))
365 if(courantmaxtots>smalldouble) dtnew=min(dtnew,
courantpar/courantmaxtots)
366 if(
fld_debug)print *,
'DT limit RADFORCE CSRAD=',dtnew
379 integer,
intent(in) :: ixi^
l, ixo^
l
380 double precision,
intent(in) :: w(ixi^s, 1:nw)
381 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
382 double precision,
intent(out) :: fld_kappa(ixi^s)
385 double precision :: rho0,temp0,n
386 double precision :: temp(ixi^s)
389 case(
'const',
'thomson')
393 {
do ix^
d=ixomin^
d,ixomax^
d\ }
401 call mpistop(
"special opacity not defined")
405 call mpistop(
"Doesn't know opacity law")
413 integer,
intent(in) :: ixi^
l, ixo^
l
414 double precision,
intent(in) :: w(ixi^s, 1:nw)
415 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
416 double precision,
intent(out):: rad_pressure(ixi^s,1:
ndim,1:
ndim)
418 double precision :: eddington_tensor(ixi^s,1:
ndim,1:
ndim)
419 double precision :: lambda(ixi^s),fld_r(ixi^s)
424 print *,
'In get_radPress with nth=',nth,
' on ixO=',ixo^
l
425 print *,
'Max and Min value of fe'
426 print *,maxval(eddington_tensor(ixo^s,1:
ndim,1:
ndim))
427 print *,minval(eddington_tensor(ixo^s,1:
ndim,1:
ndim))
428 print *,
'Max and Min value of Erad'
429 print *,maxval(w(ixo^s,iw_r_e))
430 print *,minval(w(ixo^s,iw_r_e))
431 print *,
'End get_radPress'
435 rad_pressure(ixo^s,i,j) = eddington_tensor(ixo^s,i,j)*w(ixo^s,iw_r_e)
448 integer,
intent(in) :: ixi^
l,ixo^
l,nth
449 double precision,
intent(in) :: w(ixi^s,1:nw)
450 double precision,
intent(in) :: x(ixi^s,1:
ndim)
451 double precision,
intent(out) :: fld_r(ixi^s),fld_lambda(ixi^s)
453 double precision :: wprim(ixi^s,1:nw)
455 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
469 integer,
intent(in) :: ixi^
l,ixo^
l,nth
470 double precision,
intent(in) :: w(ixi^s,1:nw)
471 double precision,
intent(in) :: x(ixi^s,1:
ndim)
472 double precision,
intent(out) :: fld_r(ixi^s),fld_lambda(ixi^s)
474 integer :: idir, ix^
d
475 double precision :: kappa(ixi^s),normgrad2(ixi^s)
476 double precision :: grad_r_e(ixi^s)
481 fld_lambda(ixo^s) = 1.d0/3.d0
485 normgrad2(ixo^s) = zero
487 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,grad_r_e,nth)
488 normgrad2(ixo^s) = normgrad2(ixo^s)+grad_r_e(ixo^s)**2
493 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
494 where(normgrad2(ixo^s)<smalldouble**2)
497 fld_lambda(ixo^s) = 1.0d0/3.0d0
499 fld_lambda(ixo^s) = one/fld_r(ixo^s)
504 normgrad2(ixo^s) = zero
506 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,grad_r_e,nth)
507 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_e(ixo^s)**2
510 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
513 fld_lambda(ixo^s) = (2.d0+fld_r(ixo^s))/(6.d0+3*fld_r(ixo^s)+fld_r(ixo^s)**2)
517 normgrad2(ixo^s) = zero
519 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,grad_r_e,nth)
520 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_e(ixo^s)**2
523 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
526 {
do ix^
d = ixomin^
d,ixomax^
d\ }
527 if(fld_r(ix^
d) .lt. 3.d0/2.d0)
then
528 fld_lambda(ix^
d) = 2.d0/(3.d0+dsqrt(9.d0+12.d0*fld_r(ix^
d)**2))
530 fld_lambda(ix^
d) = 1.d0/(1.d0+fld_r(ix^
d)+dsqrt(1.d0+2.d0*fld_r(ix^
d)))
535 call mpistop(
"special fluxlimiter not defined")
539 call mpistop(
'Fluxlimiter unknown')
551 integer,
intent(in) :: ixi^
l, ixo^
l
552 double precision,
intent(in) :: w(ixi^s, 1:nw)
553 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
554 double precision,
intent(out) :: rad_flux(ixi^s, 1:
ndim)
556 integer :: idir,nth_for_fld
557 double precision :: cc
558 double precision :: wprim(ixi^s,1:nw)
559 double precision :: grad_r_e(ixi^s)
560 double precision :: kappa(ixi^s), lambda(ixi^s), fld_r(ixi^s)
562 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
572 call gradient(wprim(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,grad_r_e,nth_for_fld)
573 rad_flux(ixo^s, idir) = -cc*lambda(ixo^s)/(kappa(ixo^s)*wprim(ixo^s,iw_rho))*grad_r_e(ixo^s)
582 integer,
intent(in) :: ixI^L, ixO^L, nth
583 double precision,
intent(in) :: w(ixI^S, 1:nw)
584 double precision,
intent(in) :: x(ixI^S, 1:ndim)
585 double precision,
intent(out) :: eddington_tensor(ixI^S,1:ndim,1:ndim)
586 double precision,
intent(out) :: lambda(ixI^S),fld_R(ixI^S)
589 double precision :: normgrad2(ixI^S)
590 double precision :: tmp(ixI^S),grad_r_e(ixI^S,1:ndim)
591 double precision :: nn_regularized(ixI^S,1:ndim,1:ndim)
593 normgrad2(ixo^s) = zero
595 call gradient(w(ixi^s, iw_r_e),ixi^l,ixo^l,idir,tmp,nth)
596 grad_r_e(ixo^s,idir)=tmp(ixo^s)
597 normgrad2(ixo^s)=normgrad2(ixo^s)+tmp(ixo^s)**2
602 nn_regularized(ixo^s,idir,jdir)=(grad_r_e(ixo^s,idir)*grad_r_e(ixo^s,jdir)+smalldouble**2)/(normgrad2(ixo^s)+smalldouble**2)
604 nn_regularized(ixo^s,idir,jdir)=(grad_r_e(ixo^s,idir)*grad_r_e(ixo^s,jdir))/(normgrad2(ixo^s)+smalldouble**2)
611 tmp(ixo^s) = lambda(ixo^s)+(lambda(ixo^s)*fld_r(ixo^s))**2
614 eddington_tensor(ixo^s,idir,idir) = half*(one-tmp(ixo^s))
619 if(idir .ne. jdir) eddington_tensor(ixo^s,idir,jdir) = zero
621 eddington_tensor(ixo^s,idir,jdir) = eddington_tensor(ixo^s,idir,jdir)+&
622 half*(3.d0*tmp(ixo^s)-one)*nn_regularized(ixo^s,idir,jdir)
634 type(state),
target :: psa(max_blocks)
635 type(state),
target :: psb(max_blocks)
636 double precision,
intent(in) :: qdt
637 double precision,
intent(in) :: qtC
638 double precision,
intent(in) :: dtfactor
640 integer,
parameter :: max_its = 100
641 integer :: n,ixO^L,ix^D
642 double precision :: res, max_residual, lambda, fac
643 integer :: iigrid, igrid
648 print *,
'skipping implicit solve: dt too small!'
655 call mg_set_methods(
mg)
656 if(.not.
mg%is_allocated)
call mpistop(
"multigrid tree not allocated yet")
663 lambda = 1.d0/(dtfactor *qdt)
665 call vhelmholtz_set_lambda(lambda)
674 call mg_restrict(
mg, mg_iveps)
675 call mg_fill_ghost_cells(
mg, mg_iveps)
677 call mg_fas_fmg(
mg, .true., max_res=res)
679 if(res < max_residual)
exit
680 call mg_fas_vcycle(
mg, max_res=res)
682 if((res .lt. 0.d0) .or. (res .ge. max_residual))
then
683 if (
mg%my_rank == 0)
then
684 write(*,*)
it,
' residual from MG ', res
685 write(*,*)
it,
' max_residual in MG ', max_residual
686 write(*,*)
it,
' qdt in MG ', qdt,
' versus dtmin=',
dtmin
687 write(*,*)
it,
' dtfactor*qdt in MG ', qdt*dtfactor
689 call mpistop(
"no convergence in MG")
691 if(res .eq. 0.d0)
then
692 if (
mg%my_rank == 0)
write(*,*)
it,
'ZERO residual from MG ', res
699 do iigrid=1,igridstail; igrid=igrids(iigrid);
700 {
do ix^d = ixomin^d,ixomax^d\ }
701 psa(igrid)%w(ix^d,iw_r_e)= max(psa(igrid)%w(ix^d,iw_r_e),
small_r_e)
710 type(state),
target :: psa(max_blocks)
712 integer :: iigrid, igrid, ixO^L
717 do iigrid=1,igridstail; igrid=igrids(iigrid);
727 type(state),
target :: psa(max_blocks)
728 double precision,
intent(in) :: qtC
729 integer :: iigrid, igrid
736 do iigrid=1,igridstail; igrid=igrids(iigrid);
747 integer,
intent(in) :: ixI^L, ixO^L
748 double precision,
intent(inout) :: w(ixI^S, 1:nw)
750 double precision :: divF(ixI^S)
751 integer :: idir, jxO^L, hxO^L
753 if(.not.
slab)
call mpistop(
"laplacian coded up for uniform cartesian grid")
760 hxo^l=ixo^l-
kr(idir,^
d);
761 jxo^l=ixo^l+
kr(idir,^
d);
762 divf(ixo^s) = divf(ixo^s) + &
774 w(ixo^s,iw_r_e) = divf(ixo^s)
784 integer,
intent(in) :: ixi^
l, ixo^
l
785 double precision,
intent(in) :: x(ixi^s,1:
ndim)
786 double precision,
intent(inout) :: w(ixi^s,1:nw)
788 double precision :: cc
789 double precision :: wprim(ixi^s,1:nw)
790 double precision :: kappa(ixi^s),lambda(ixi^s),fld_r(ixi^s)
793 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
799 w(ixo^s,
i_diff_mg) = cc*lambda(ixo^s)/(kappa(ixo^s)*wprim(ixo^s,iw_rho))
801 if(minval(w(ixo^s,
i_diff_mg))<smalldouble)
then
802 print *,
'min diffcoef=',minval(w(ixo^s,
i_diff_mg))
803 call mpistop(
"too small diffusion coefficient")
805 if(maxval(w(ixo^s,
i_diff_mg))>bigdouble)
call mpistop(
"too large diffusion coefficient")
808 print *,
'setting diffcoefs with data on',ixi^
l
809 print *,
'min diffcoef=',minval(w(ixo^s,
i_diff_mg))
810 print *,
'min lambda kappa rho fld_R'
811 print *,minval(lambda(ixo^s))
812 print *,minval(kappa(ixo^s))
813 print *,minval(wprim(ixo^s,iw_rho))
814 print *,minval(fld_r(ixo^s))
815 print *,
'max diffcoef=',maxval(w(ixo^s,
i_diff_mg))
816 print *,
'max lambda kappa rho fld_R'
817 print *,maxval(lambda(ixo^s))
818 print *,maxval(kappa(ixo^s))
819 print *,maxval(wprim(ixo^s,iw_rho))
820 print *,maxval(fld_r(ixo^s))
821 print *,
'done setting diffcoefs in slot',
i_diff_mg,
' on range',ixo^
l
830 double precision,
intent(in) :: c0, c1
831 double precision,
intent(in) :: E_rad
832 double precision,
intent(inout) :: e_gas
833 double precision :: bisect_a, bisect_b, bisect_c
834 integer :: n, max_its
839 bisect_b = min(dabs(c0/c1),dabs(c0)**(1.d0/4.d0))+smalldouble
840 do while(dabs(bisect_b-bisect_a) .ge.
fld_bisect_tol*min(e_gas,e_rad))
841 bisect_c = (bisect_a + bisect_b)/two
843 if(n .gt. max_its)
then
844 call mpistop(
'No convergece in bisection scheme')
867 call mpistop(
"Problem with fld bisection method")
875 if(
fld_debug)print*,
"IGNORING GAS-RAD ENERGY EXCHANGE ", c0, c1
877 call mpistop(
'issues in bisection scheme')
886 2435 e_gas = (bisect_a + bisect_b)/two
892 double precision,
intent(in) :: c0, c1
893 double precision,
intent(in) :: E_rad
894 double precision,
intent(inout) :: e_gas
895 double precision :: xval, yval, der, deltax
911 if(
fld_debug)print*,
'skip to bisection algorithm'
922 double precision,
intent(in) :: c0, c1
923 double precision,
intent(in) :: E_rad
924 double precision,
intent(inout) :: e_gas
925 double precision :: xval, yval, der, dder, deltax
939 deltax = -two*yval*der/(two*der**2 - yval*dder)
943 if(
fld_debug)print*,
'skip to Newton algorithm'
955 double precision,
intent(in) :: e_gas
956 double precision,
intent(in) :: c0, c1
957 double precision :: val
959 val = e_gas**4.d0 + c1*e_gas - c0
965 double precision,
intent(in) :: e_gas
966 double precision,
intent(in) :: c0, c1
967 double precision :: der
969 der = 4.d0*e_gas**3.d0 + c1
975 double precision,
intent(in) :: e_gas
976 double precision,
intent(in) :: c0, c1
977 double precision :: dder
979 dder = 4.d0*3.d0*e_gas**2.d0
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter half
double precision, parameter one
double precision, parameter const_c
double precision, parameter const_rad_a
Module for flux limited diffusion (FLD)-approximation in Radiation-(Magneto)hydrodynamics simulations...
subroutine, public fld_get_radpress(w, x, ixil, ixol, rad_pressure)
Returns Radiation Pressure as tensor NOTE: w is primitive on entry.
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 radiative Energy diffusion.
subroutine fld_params_read(files)
public methods these are called in mod_hd_phys or mod_mhd_phys
subroutine update_diffcoeff(psa)
character(len=8) fld_opacity_law
switches for opacity
subroutine, public fld_get_diffcoef_central(w, x, ixil, ixol)
Calculates cell-centered diffusion coefficient to be used in multigrid.
subroutine evaluate_diffterm_onegrid(ixil, ixol, w)
inplace update of psa==>F_im(psa)
double precision function ddpolynomial_bisection(e_gas, c0, c1)
Evaluate second derivative of polynomial at argument e_gas.
character(len=16) fld_fluxlimiter
flux limiter choice
subroutine, public fld_get_radflux(w, x, ixil, ixol, rad_flux)
Calculate Radiation Flux (use of cgs units) NOTE: only for diagnostics purposes (w conservative on en...
subroutine bisection_method(e_gas, e_rad, c0, c1)
Find the root of the 4th degree polynomial using the bisection method.
character(len=50) fld_opal_table
double precision, public fld_kappa0
Opacity per unit of unit_density.
subroutine, public fld_get_fluxlimiter(w, x, ixil, ixol, fld_lambda, fld_r, nth)
This subroutine calculates flux limiter lambda according to fld_fluxlimiter It also calculates fld_R ...
double precision function polynomial_bisection(e_gas, c0, c1)
Evaluate polynomial at argument e_gas.
subroutine fld_get_eddington(w, x, ixil, ixol, eddington_tensor, lambda, fld_r, nth)
Calculate Eddington-tensor (where w is primitive) also feeds back the flux limiter lambda and R.
subroutine, public fld_get_opacity_prim(w, x, ixil, ixol, fld_kappa)
Sets the opacity in the w-array by calling mod_opal_opacity NOTE: assumes primitives in w NOTE: assum...
subroutine, public fld_init(he_abundance, r_gamma)
Initialising FLD-module Read opacities Initialise Multigrid and adimensionalise kappa.
subroutine fld_set_mg_bounds
Set the boundaries for the diffusion of E.
subroutine halley_method(e_gas, e_rad, c0, c1)
Find the root of the 4th degree polynomial using the Halley method.
subroutine fld_evaluate_implicit(qtc, psa)
inplace update of psa==>F_im(psa)
subroutine, public add_fld_rad_force(qdt, ixil, ixol, wct, wctprim, w, x, 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...
character(len=8) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
subroutine newton_method(e_gas, e_rad, c0, c1)
Find the root of the 4th degree polynomial using the Newton method.
logical fld_debug
switches for local debug purposes
logical fld_radforce_split
source split for energy interact and radforce:
double precision function dpolynomial_bisection(e_gas, c0, c1)
Evaluate first derivative of polynomial at argument e_gas.
subroutine, public fld_get_fluxlimiter_prim(w, x, ixil, ixol, fld_lambda, fld_r, nth)
This subroutine calculates flux limiter lambda according to fld_fluxlimiter It also calculates fld_R ...
subroutine, public fld_radforce_get_dt(w, ixil, ixol, dtnew, dxd, x)
get dt limit for radiation force and FLD explicit source additions NOTE: w is primitive on entry
integer i_diff_mg
diffusion coefficient for multigrid method
subroutine fld_implicit_update(dtfactor, qdt, qtc, psa, psb)
Calling all subroutines to perform the multigrid method Communicates rad_e and diff_coeff to multigri...
integer nth_for_diff_mg
diffusion coefficient stencil control
Module with basic grid data structures.
Module with geometry-related routines (e.g., divergence, curl)
subroutine gradient(q, ixil, ixol, idir, gradq, nth_in)
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
integer, parameter bc_noinflow
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
integer, parameter bc_asymm
double precision global_time
The global simulation time.
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer it
Number of time steps taken.
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
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.
double precision, dimension(:), allocatable, parameter d
double precision courantpar
The Courant (CFL) number used for the simulation.
integer ixm
the mesh range of a physical block without ghost cells
logical slab
Cartesian geometry or not.
integer, parameter bc_periodic
integer, parameter bc_special
boundary condition types
double precision unit_velocity
Physical scaling factor for velocity.
logical time_advance
do time evolving
double precision small_r_e
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision unit_temperature
Physical scaling factor for temperature.
integer, parameter bc_cont
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
integer, parameter bc_symm
logical fix_small_values
fix small values with average or replace methods
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical use_multigrid
Use multigrid (only available in 2D and 3D)
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
double precision dtmin
Stop the simulation when the time step becomes smaller than this value.
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(tabledir, set_custom_tabledir)
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_convert), pointer phys_to_primitive
procedure(sub_get_tgas), pointer phys_get_tgas
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
procedure(sub_get_rfactor), pointer phys_get_rfactor
integer phys_wider_stencil
To use wider stencils in flux calculations. A value of 1 will extend it by one cell in both direction...
procedure(sub_set_mg_bounds), pointer phys_set_mg_bounds
procedure(sub_implicit_update), pointer phys_implicit_update
procedure(sub_get_csrad2), pointer phys_get_csrad2
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_mg_bc), pointer usr_special_mg_bc
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...