61 character(len=*),
intent(in) :: files(:)
69 open(
unitpar, file=trim(files(n)), status=
"old")
70 read(
unitpar, fld_list,
end=111)
85 double precision,
intent(in) :: r_gamma
96 call mpistop(
"please set the constant opacity to a reasonable value")
100 call mpistop(
"convergence tolerance for root solver too strict")
105 call mpistop(
"convergence tolerance for MG solver too strict")
114 call mpistop(
"nth_for_diff_mg must be 1 or 2")
125 mg%operator_type = mg_vhelmholtz
127 mg%smoother_type = mg_smoother_gsrb
133 if(
si_unit)
call mpistop(
"adjust opal module with SI-cgs conversions for SI - or use cgs!")
151 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
152 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
155 mg%bc(ib, mg_iphi)%bc_type = mg_bc_dirichlet
156 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
159 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
160 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
165 print *,
'Special boundary for Erad needs specific user-set MG BC treatment'
166 print *,
' and this could be through usr_special_mg_bc call'
172 call mpistop(
"divE_multigrid warning: unknown b.c. ")
176 mg%bc(ib, mg_iveps)%bc_type = mg_bc_neumann
177 mg%bc(ib, mg_iveps)%bc_value = 0.0_dp
191 integer,
intent(in) :: ixi^
l, ixo^
l
192 double precision,
intent(in) :: qdt, x(ixi^s,1:
ndim)
193 double precision,
intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw)
194 double precision,
intent(inout) :: w(ixi^s,1:nw)
195 logical,
intent(in) :: qsourcesplit
196 logical,
intent(inout) :: active
198 integer :: idir,jdir,nth_for_fld,ix^
d
199 double precision,
dimension(ixI^S) :: a1,a2,a3,c0,c1,kappa
200 double precision,
dimension(ixI^S) :: e_gas,e_rad,tmp
201 double precision,
dimension(ixI^S,1:ndim,1:ndim) :: div_v,edd
210 call gradient(wctprim(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,tmp,nth_for_fld)
213 tmp(ixo^s) = -a1(ixo^s)*tmp(ixo^s)
215 w(ixo^s,iw_mom(idir)) = w(ixo^s,iw_mom(idir))+ qdt*tmp(ixo^s)
217 w(ixo^s,iw_e) = w(ixo^s,iw_e) + qdt*wctprim(ixo^s,iw_mom(idir))*tmp(ixo^s)
225 call gradient(wctprim(ixi^s,iw_mom(jdir)),ixi^
l,ixo^
l,idir,tmp)
226 div_v(ixo^s,idir,jdir) = tmp(ixo^s)
231 a3(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1)
234 a3(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1) &
235 + div_v(ixo^s,1,2)*edd(ixo^s,1,2) &
236 + div_v(ixo^s,2,1)*edd(ixo^s,2,1) &
237 + div_v(ixo^s,2,2)*edd(ixo^s,2,2)
240 a3(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1) &
241 + div_v(ixo^s,1,2)*edd(ixo^s,1,2) &
242 + div_v(ixo^s,1,3)*edd(ixo^s,1,3) &
243 + div_v(ixo^s,2,1)*edd(ixo^s,2,1) &
244 + div_v(ixo^s,2,2)*edd(ixo^s,2,2) &
245 + div_v(ixo^s,2,3)*edd(ixo^s,2,3) &
246 + div_v(ixo^s,3,1)*edd(ixo^s,3,1) &
247 + div_v(ixo^s,3,2)*edd(ixo^s,3,2) &
248 + div_v(ixo^s,3,3)*edd(ixo^s,3,3)
253 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)
254 if(
allocated(iw_mag))
then
255 e_gas(ixo^s) = e_gas(ixo^s)-
half*sum(w(ixo^s,iw_mag(:))**2,dim=
ndim+1)
257 e_rad(ixo^s) = w(ixo^s,iw_r_e)
260 {
do ix^db= ixomin^db,ixomax^db\}
262 write(*,*)
"Error in FLD add_fld_rad_force: small value"
263 write(*,*)
"of internal or radiation energy density before exchange"
265 write(*,*)
"Location: ", x(ix^
d,:)
266 write(*,*)
"Cell number: ", ix^
d
267 write(*,*)
"internal energy density is=",e_gas(ix^
d),
" versus small_e=",
small_e
268 write(*,*)
"radiation energy density is=",e_rad(ix^
d),
" versus small_r_e=",
small_r_e
269 call mpistop(
"FLD error:May need to turn on fixes")
274 if(fix_small_values)
then
275 {
do ix^d = ixomin^d,ixomax^d\ }
276 e_gas(ix^d) = max(e_gas(ix^d),small_e)
277 e_rad(ix^d) = max(e_rad(ix^d),small_r_e)
283 call phys_get_rfactor(wct,x,ixi^l,ixo^l,tmp)
284 a1(ixo^s) = qdt*kappa(ixo^s)*c_norm*arad_norm*(
fld_gamma-one)**4/(wct(ixo^s,iw_rho)**3*tmp(ixo^s)**4)
285 a2(ixo^s) = c_norm*kappa(ixo^s)*wct(ixo^s,iw_rho)*qdt
286 a3(ixo^s) = a3(ixo^s)*qdt
288 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))
289 c1(ixo^s) = (one+a2(ixo^s)+a3(ixo^s))/a1(ixo^s)/(one+a3(ixo^s))
292 {
do ix^d = ixomin^d,ixomax^d\}
301 call mpistop(
'root-method not known')
305 e_rad(ixo^s) = (a1(ixo^s)*e_gas(ixo^s)**4.d0+e_rad(ixo^s))/(one+a2(ixo^s)+a3(ixo^s))
307 if(check_small_values.and..not.fix_small_values)
then
308 {
do ix^db= ixomin^db,ixomax^db\}
309 if(e_gas(ix^d)<small_e.or.e_rad(ix^d)<small_r_e)
then
310 write(*,*)
"Error in FLD add_fld_rad_force: small value"
311 write(*,*)
"of internal or radiation energy density after exchange"
312 write(*,*)
"Iteration: ", it,
" Time: ", global_time
313 write(*,*)
"Location: ", x(ix^d,:)
314 write(*,*)
"Cell number: ", ix^d
315 write(*,*)
"internal energy density is=",e_gas(ix^d),
" versus small_e=",small_e
316 write(*,*)
"radiation energy density is=",e_rad(ix^d),
" versus small_r_e=",small_r_e
317 call mpistop(
"FLD error:May need to turn on fixes")
322 if(fix_small_values)
then
323 {
do ix^d = ixomin^d,ixomax^d\ }
324 e_gas(ix^d) = max(e_gas(ix^d),small_e)
325 e_rad(ix^d) = max(e_rad(ix^d),small_r_e)
330 w(ixo^s,iw_e) = e_gas(ixo^s)
331 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)
332 if(
allocated(iw_mag))
then
333 w(ixo^s,iw_e) = w(ixo^s,iw_e)+half*sum(w(ixo^s,iw_mag(:))**2,dim=ndim+1)
336 w(ixo^s,iw_r_e) = e_rad(ixo^s)
348 integer,
intent(in) :: ixi^
l, ixo^
l
349 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim), w(ixi^s,1:nw)
350 double precision,
intent(inout) :: dtnew
352 integer :: idim,idims,nth_for_fld
353 double precision :: dxinv(1:
ndim), max_grav
354 double precision :: lambda(ixi^s),fld_r(ixi^s)
355 double precision :: tmp(ixi^s)
356 double precision :: cmax(ixi^s),cmaxtot(ixi^s),courantmaxtots
358 if(
fld_debug)print *,
'DT limit on entry to radforce_get_dt=',dtnew
362 ^
d&dxinv(^
d)=one/
dx^
d;
364 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idim,tmp,nth_for_fld)
365 max_grav = maxval(dabs(-lambda(ixo^s)*tmp(ixo^s)/w(ixo^s,iw_rho)))
366 max_grav = max(max_grav, epsilon(1.0d0))
367 dtnew = min(dtnew, 1.0d0 / dsqrt(max_grav * dxinv(idim)))
371 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idim,tmp,nth_for_fld)
372 max_grav = maxval(dabs(-lambda(ixo^s)*tmp(ixo^s)/w(ixo^s,iw_rho))/
block%ds(ixo^s,idim))
373 max_grav = max(max_grav, epsilon(1.0d0))
374 dtnew = min(dtnew, 1.0d0 / dsqrt(max_grav))
378 if(
fld_debug)print *,
'DT limit after RADFORCE eff grav=',dtnew
383 ^
d&dxinv(^
d)=one/
dx^
d;
385 cmax(ixo^s)=dabs(w(ixo^s,iw_mom(idims)))+dsqrt(tmp(ixo^s))
387 cmaxtot(ixo^s)=cmax(ixo^s)*dxinv(idims)
389 cmaxtot(ixo^s)=cmaxtot(ixo^s)+cmax(ixo^s)*dxinv(idims)
394 cmax(ixo^s)=dabs(w(ixo^s,iw_mom(idims)))+dsqrt(tmp(ixo^s))
396 cmaxtot(ixo^s)=cmax(ixo^s)/
block%ds(ixo^s,idims)
398 cmaxtot(ixo^s)=cmaxtot(ixo^s)+cmax(ixo^s)/
block%ds(ixo^s,idims)
403 courantmaxtots=maxval(cmaxtot(ixo^s))
404 if(courantmaxtots>smalldouble) dtnew=min(dtnew,
courantpar/courantmaxtots)
405 if(
fld_debug)print *,
'DT limit RADFORCE CSRAD=',dtnew
418 integer,
intent(in) :: ixi^
l, ixo^
l
419 double precision,
intent(in) :: w(ixi^s, 1:nw)
420 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
421 double precision,
intent(out) :: fld_kappa(ixi^s)
424 double precision :: rho0,temp0,kapp0
425 double precision :: temp(ixi^s)
434 {
do ix^
d=ixomin^
d,ixomax^
d\ }
442 call mpistop(
"special opacity not defined")
446 call mpistop(
"Doesn't know opacity law")
454 integer,
intent(in) :: ixi^
l, ixo^
l
455 double precision,
intent(in) :: w(ixi^s, 1:nw)
456 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
457 double precision,
intent(out):: rad_pressure(ixi^s,1:
ndim,1:
ndim)
459 double precision :: eddington_tensor(ixi^s,1:
ndim,1:
ndim)
460 double precision :: lambda(ixi^s),fld_r(ixi^s)
466 print *,
'In get_radPress with nth=',nth,
' on ixO=',ixo^
l
467 print *,
'Max and Min value of fe'
468 print *,maxval(eddington_tensor(ixo^s,1:
ndim,1:
ndim))
469 print *,minval(eddington_tensor(ixo^s,1:
ndim,1:
ndim))
470 print *,
'Max and Min value of Erad'
471 print *,maxval(w(ixo^s,iw_r_e))
472 print *,minval(w(ixo^s,iw_r_e))
473 print *,
'End get_radPress'
477 rad_pressure(ixo^s,i,j) = eddington_tensor(ixo^s,i,j)*w(ixo^s,iw_r_e)
490 integer,
intent(in) :: ixi^
l,ixo^
l,nth
491 double precision,
intent(in) :: w(ixi^s,1:nw)
492 double precision,
intent(in) :: x(ixi^s,1:
ndim)
493 double precision,
intent(out) :: fld_r(ixi^s),fld_lambda(ixi^s)
495 double precision :: wprim(ixi^s,1:nw)
497 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
511 integer,
intent(in) :: ixi^
l,ixo^
l,nth
512 double precision,
intent(in) :: w(ixi^s,1:nw)
513 double precision,
intent(in) :: x(ixi^s,1:
ndim)
514 double precision,
intent(out) :: fld_r(ixi^s),fld_lambda(ixi^s)
516 integer :: idir, ix^
d
517 double precision :: kappa(ixi^s),normgrad2(ixi^s)
518 double precision :: grad_r_e(ixi^s)
523 fld_lambda(ixo^s) = 1.d0/3.d0
527 normgrad2(ixo^s) = zero
529 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,grad_r_e,nth)
530 normgrad2(ixo^s) = normgrad2(ixo^s)+grad_r_e(ixo^s)**2
535 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
536 where(normgrad2(ixo^s)<smalldouble**2)
539 fld_lambda(ixo^s) = 1.0d0/3.0d0
541 fld_lambda(ixo^s) = one/fld_r(ixo^s)
546 normgrad2(ixo^s) = zero
548 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,grad_r_e,nth)
549 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_e(ixo^s)**2
552 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
555 fld_lambda(ixo^s) = (2.d0+fld_r(ixo^s))/(6.d0+3*fld_r(ixo^s)+fld_r(ixo^s)**2)
559 normgrad2(ixo^s) = zero
561 call gradient(w(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,grad_r_e,nth)
562 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_e(ixo^s)**2
565 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
568 {
do ix^
d = ixomin^
d,ixomax^
d\ }
569 if(fld_r(ix^
d) .lt. 3.d0/2.d0)
then
570 fld_lambda(ix^
d) = 2.d0/(3.d0+dsqrt(9.d0+12.d0*fld_r(ix^
d)**2))
572 fld_lambda(ix^
d) = 1.d0/(1.d0+fld_r(ix^
d)+dsqrt(1.d0+2.d0*fld_r(ix^
d)))
577 call mpistop(
"special fluxlimiter not defined")
581 call mpistop(
'Fluxlimiter unknown')
593 integer,
intent(in) :: ixi^
l, ixo^
l
594 double precision,
intent(in) :: w(ixi^s, 1:nw)
595 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
596 double precision,
intent(out) :: rad_flux(ixi^s, 1:
ndim)
598 integer :: idir,nth_for_fld
599 double precision :: wprim(ixi^s,1:nw)
600 double precision :: grad_r_e(ixi^s)
601 double precision :: kappa(ixi^s), lambda(ixi^s), fld_r(ixi^s)
603 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
613 call gradient(wprim(ixi^s,iw_r_e),ixi^
l,ixo^
l,idir,grad_r_e,nth_for_fld)
614 rad_flux(ixo^s,idir)=-(
c_norm*lambda(ixo^s)/(kappa(ixo^s)*wprim(ixo^s,iw_rho)))*grad_r_e(ixo^s)
623 integer,
intent(in) :: ixI^L, ixO^L, nth
624 double precision,
intent(in) :: w(ixI^S, 1:nw)
625 double precision,
intent(in) :: x(ixI^S, 1:ndim)
626 double precision,
intent(out) :: eddington_tensor(ixI^S,1:ndim,1:ndim)
627 double precision,
intent(out) :: lambda(ixI^S),fld_R(ixI^S)
630 double precision :: normgrad2(ixI^S)
631 double precision :: tmp(ixI^S),grad_r_e(ixI^S,1:ndim)
632 double precision :: nn_regularized(ixI^S,1:ndim,1:ndim)
634 normgrad2(ixo^s) = zero
636 call gradient(w(ixi^s, iw_r_e),ixi^l,ixo^l,idir,tmp,nth)
637 grad_r_e(ixo^s,idir)=tmp(ixo^s)
638 normgrad2(ixo^s)=normgrad2(ixo^s)+tmp(ixo^s)**2
643 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)
645 nn_regularized(ixo^s,idir,jdir)=(grad_r_e(ixo^s,idir)*grad_r_e(ixo^s,jdir))/(normgrad2(ixo^s)+smalldouble**2)
652 tmp(ixo^s) = lambda(ixo^s)+(lambda(ixo^s)*fld_r(ixo^s))**2
655 eddington_tensor(ixo^s,idir,idir) = half*(one-tmp(ixo^s))
660 if(idir .ne. jdir) eddington_tensor(ixo^s,idir,jdir) = zero
662 eddington_tensor(ixo^s,idir,jdir) = eddington_tensor(ixo^s,idir,jdir)+&
663 half*(3.d0*tmp(ixo^s)-one)*nn_regularized(ixo^s,idir,jdir)
677 type(state),
target :: psa(max_blocks)
678 type(state),
target :: psb(max_blocks)
679 double precision,
intent(in) :: qdt
680 double precision,
intent(in) :: qtC
681 double precision,
intent(in) :: dtfactor
683 integer,
parameter :: max_its = 100
684 integer :: n,ixO^L,ix^D
685 double precision :: res, max_residual, mg_lambda, fac
686 double precision :: wmax(nw),wmin(nw)
687 double precision :: wmaxb(nw),wminb(nw)
688 integer :: iigrid, igrid
694 print *,
'skipping implicit solve: dt too small!'
704 mg_lambda = 1.d0/(dtfactor *qdt)
715 print *,
'at start of MG solver, we have fld_diff_tol =',
fld_diff_tol
716 print *,
'at start of MG solver, we have LHS E_rad range as :',wmax(iw_r_e),wmin(iw_r_e)
717 print *,
'at start of MG solver, we have Diff coeff range as:',wmax(
i_diff_mg),wmin(
i_diff_mg)
718 print *,
'at start of MG solver, we have density range as :',wmax(iw_rho),wmin(iw_rho)
719 print *,
'at start of MG solver, we have RHS E_rad range as :',wmaxb(iw_r_e),wminb(iw_r_e)
720 print *,
'at start of MG solver, we have qdt as',qdt,
' and max_residual=',max_residual
721 print *,
'at start of MG solver, we have dtfactor as',dtfactor,
' or 1/dt (or mg_lambda)=',mg_lambda
722 print *,
'at start of MG solver, ratio coeffs on level 1 =',wmax(
i_diff_mg)/(mg_lambda*
dx(1,1)**2)
727 call mg_set_methods(
mg)
728 if(.not.
mg%is_allocated)
call mpistop(
"multigrid tree not allocated yet")
735 call vhelmholtz_set_lambda(mg_lambda)
742 call mg_copy_to_tree(iw_r_e, mg_irhs, factor=-mg_lambda, state_from=psb)
745 call mg_restrict(
mg, mg_iveps)
746 call mg_fill_ghost_cells(
mg, mg_iveps)
748 call mg_fas_fmg(
mg, .true., max_res=res)
749 if(
fld_debug.and.
mype==0)print *,
'MG residual obtained with FMG is =',res
751 if(res < max_residual)
exit
752 call mg_fas_vcycle(
mg, max_res=res)
753 if(
fld_debug.and.
mype==0)print *,
'MG residual obtained with Vcycle =',res,
' at step =',n
755 if(
fld_debug.and.
mype==0)print *,
'FINAL MG residual obtained is =',res
756 if(res .ge. max_residual)
then
758 write(*,*)
it,
' residual from MG ', res
759 write(*,*)
it,
' max_residual in MG ', max_residual
760 write(*,*)
it,
' qdt in MG ', qdt,
' versus dtmin=',
dtmin
761 write(*,*)
it,
' dtfactor*qdt in MG ', qdt*dtfactor
764 call mpistop(
"no convergence in MG")
766 if(
mype==0)
write(*,*)
'WARNING for it=',
it,
' NO CONVERGENCE IN MG but still carry on'
774 do iigrid=1,igridstail; igrid=igrids(iigrid);
775 {
do ix^db= ixomin^db,ixomax^db\}
776 if(psa(igrid)%w(ix^d,iw_r_e)<
small_r_e)
then
777 write(*,*)
"Error in FLD fld_implicit_update: small value"
778 write(*,*)
"of radiation energy density after MG"
780 write(*,*)
"Location: ", psa(igrid)%x(ix^d,:),
" on grid",igrid
781 write(*,*)
"Cell number: ", ix^d
782 write(*,*)
"radiation energy density is=",psa(igrid)%w(ix^d,iw_r_e),
" versus small_r_e=",
small_r_e
783 call mpistop(
"FLD error:May need to turn on fixes")
789 if(fix_small_values)
then
791 do iigrid=1,igridstail; igrid=igrids(iigrid);
792 {
do ix^d = ixomin^d,ixomax^d\ }
793 psa(igrid)%w(ix^d,iw_r_e)= max(psa(igrid)%w(ix^d,iw_r_e),small_r_e)
802 type(state),
target :: psa(max_blocks)
804 integer :: iigrid, igrid, ixO^L
809 do iigrid=1,igridstail; igrid=igrids(iigrid);
819 type(state),
target :: psa(max_blocks)
820 double precision,
intent(in) :: qtC
821 integer :: iigrid, igrid
828 do iigrid=1,igridstail; igrid=igrids(iigrid);
839 integer,
intent(in) :: ixI^L, ixO^L
840 double precision,
intent(inout) :: w(ixI^S, 1:nw)
842 double precision :: divF(ixI^S)
843 integer :: idir, jxO^L, hxO^L
845 if(.not.
slab)
call mpistop(
"laplacian coded up for uniform cartesian grid")
852 hxo^l=ixo^l-
kr(idir,^
d);
853 jxo^l=ixo^l+
kr(idir,^
d);
854 divf(ixo^s) = divf(ixo^s) + &
868 w(ixo^s,iw_r_e) = divf(ixo^s)
878 integer,
intent(in) :: ixi^
l, ixo^
l
879 double precision,
intent(in) :: x(ixi^s,1:
ndim)
880 double precision,
intent(inout) :: w(ixi^s,1:nw)
882 double precision :: wprim(ixi^s,1:nw)
883 double precision :: kappa(ixi^s),lambda(ixi^s),fld_r(ixi^s)
885 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
891 w(ixo^s,
i_diff_mg) =
c_norm*lambda(ixo^s)/(kappa(ixo^s)*wprim(ixo^s,iw_rho))
893 if(minval(w(ixo^s,
i_diff_mg))<smalldouble)
then
894 print *,
'min diffcoef=',minval(w(ixo^s,
i_diff_mg))
895 call mpistop(
"too small diffusion coefficient")
897 if(maxval(w(ixo^s,
i_diff_mg))>bigdouble)
call mpistop(
"too large diffusion coefficient")
900 print *,
'setting diffcoefs with data on',ixi^
l
901 print *,
'min diffcoef=',minval(w(ixo^s,
i_diff_mg))
902 print *,
'min lambda kappa rho fld_R'
903 print *,minval(lambda(ixo^s))
904 print *,minval(kappa(ixo^s))
905 print *,minval(wprim(ixo^s,iw_rho))
906 print *,minval(fld_r(ixo^s))
907 print *,
'max diffcoef=',maxval(w(ixo^s,
i_diff_mg))
908 print *,
'max lambda kappa rho fld_R'
909 print *,maxval(lambda(ixo^s))
910 print *,maxval(kappa(ixo^s))
911 print *,maxval(wprim(ixo^s,iw_rho))
912 print *,maxval(fld_r(ixo^s))
913 print *,
'done setting diffcoefs in slot',
i_diff_mg,
' on range',ixo^
l
922 double precision,
intent(in) :: c0, c1
923 double precision,
intent(inout) :: e_gas
924 double precision :: bisect_a, bisect_b, bisect_c
925 integer :: n, max_its
930 bisect_b = min(dabs(c0/c1),dabs(c0)**(1.d0/4.d0))+smalldouble
932 bisect_c = (bisect_a + bisect_b)/two
934 if(n .gt. max_its)
then
935 call mpistop(
'No convergece in bisection scheme')
958 call mpistop(
"Problem with fld bisection method")
966 if(
fld_debug)print*,
"IGNORING GAS-RAD ENERGY EXCHANGE ", c0, c1
968 call mpistop(
'issues in bisection scheme')
977 2435 e_gas = (bisect_a + bisect_b)/two
983 double precision,
intent(in) :: c0, c1
984 double precision,
intent(inout) :: e_gas
985 double precision :: xval, yval, der, deltax
1000 if(ii .gt. 1d3)
then
1001 if(
fld_debug)print*,
'skip to bisection algorithm'
1012 double precision,
intent(in) :: c0, c1
1013 double precision,
intent(inout) :: e_gas
1014 double precision :: xval, yval, der, dder, deltax
1028 deltax = -two*yval*der/(two*der**2 - yval*dder)
1029 xval = xval + deltax
1031 if(ii .gt. 1d3)
then
1032 if(
fld_debug)print*,
'skip to Newton algorithm'
1043 double precision,
intent(in) :: e_gas
1044 double precision,
intent(in) :: c0, c1
1045 double precision :: val
1047 val = e_gas**4.d0 + c1*e_gas - c0
1053 double precision,
intent(in) :: e_gas
1054 double precision,
intent(in) :: c0, c1
1055 double precision :: der
1057 der = 4.d0*e_gas**3.d0 + c1
1063 double precision,
intent(in) :: e_gas
1064 double precision,
intent(in) :: c0, c1
1065 double precision :: dder
1067 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
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...
logical fld_force_mg_converged
switches for local debug purposes
double precision, public fld_diff_tol
Tolerance for radiative Energy diffusion.
double precision, public fld_gamma
A copy of (m)hd_gamma.
subroutine fld_params_read(files)
public methods these are called in mod_hd_phys or mod_mhd_phys
character(len=40) fld_fluxlimiter
flux limiter choice
character(len=40) fld_opal_table
subroutine update_diffcoeff(psa)
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.
subroutine, public fld_get_radflux(w, x, ixil, ixol, rad_flux)
Calculate Radiation Flux NOTE: only for diagnostics purposes (w conservative on entry) This returns c...
double precision, public fld_kappa0
Opacity value when using constant opacity.
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 ...
subroutine newton_method(e_gas, c0, c1)
Find the root of the 4th degree polynomial using the Newton method.
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...
character(len=40) fld_opacity_law
switches for opacity
character(len=40) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
subroutine fld_set_mg_bounds
Set the boundaries for the diffusion of E.
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...
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_init(r_gamma)
Initialising FLD-module Read opacities Initialise Multigrid and adimensionalise kappa.
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 bisection_method(e_gas, c0, c1)
Find the root of the 4th degree polynomial using the bisection 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...
subroutine halley_method(e_gas, c0, c1)
Find the root of the 4th degree polynomial using the Halley method.
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.
double precision const_kappae
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.
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 small_r_e
double precision c_norm
Normalised speed of light.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision unit_temperature
Physical scaling factor for temperature.
integer, parameter bc_cont
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
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.
integer refine_max_level
Maximal number of AMR levels.
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
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...