MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_fld.t
Go to the documentation of this file.
1!> Module for flux limited diffusion (FLD)-approximation in Radiation-(Magneto)hydrodynamics simulations
2!>
3!> Full description of RHD-FLD in
4!> Moens N., Sundqvist J.O., El Mellah I., Poniatowski L., Teunissen J. & Keppens R. 2022, A&A 657, A81
5!> Radiation-hydrodynamics with MPI-AMRVAC . Flux-limited diffusion
6!> doi:10.1051/0004-6361/202141023
7!>
8!> Full description for RMHD-FLD in
9!> N. Narechania, R. Keppens, A. ud-Doula, N. Moens & J. Sundqvist 2025, A&A 696, A131
10!> doi:10.1051/0004-6361/202452208
11!> Radiation-magnetohydrodynamics with MPI-AMRVAC using flux-limited diffusion
12
13module mod_fld
14 use mod_comm_lib, only: mpistop
15 use mod_geometry
16 implicit none
17 !> source split for energy interact and radforce:
18 logical :: fld_radforce_split = .false.
19 !> switch to handle photon tiring explicit or implicit
20 logical :: fld_tiring_explicit = .false.
21 !> Opacity value when using constant opacity
22 double precision, public :: fld_kappa0 = 0.0d0
23 !> Tolerance for bisection method for Energy sourceterms
24 !> This is a percentage of the minimum of gas- and radiation energy
25 double precision, public :: fld_bisect_tol = 1.d-4
26 !> Tolerance for radiative Energy diffusion
27 double precision, public :: fld_diff_tol = 1.d-4
28 !> handling ramp-up phase with explicit diffusion dt limit, slowly boosted
29 logical :: fld_slowsteps = .false.
30 double precision, public :: fld_boost_dt = 0.0d0
31 !> switches for local debug purposes
32 logical :: fld_force_mg_converged = .true.
34 !> switches for opacity
35 character(len=40) :: fld_opacity_law = 'const'
36 character(len=40) :: fld_opal_table = 'Y09800'
37 !> flux limiter choice
38 character(len=40) :: fld_fluxlimiter = 'Pomraning'
39 !> diffusion coefficient for multigrid method
40 integer :: i_diff_mg
41 !> diffusion coefficient stencil control
42 integer :: nth_for_diff_mg
43 !> Which method to find the root for the energy interaction polynomial
44 character(len=40) :: fld_interaction_method = 'Halley'
45 !> A copy of (m)hd_gamma
46 double precision, public :: fld_gamma
47
48 !> public methods
49 !> these are called in mod_hd_phys or mod_mhd_phys
50 public :: fld_init
51 public :: fld_get_radpress
53 public :: add_fld_rad_force
54 public :: fld_radforce_get_dt
55 ! these are made public for mod_usr purposes and diagnostics
56 public :: fld_get_radflux
57 public :: fld_get_fluxlimiter
59 public :: fld_get_opacity_prim
60 contains
61
62 !> Reading in fld-list parameters from .par file
63 subroutine fld_params_read(files)
66 character(len=*), intent(in) :: files(:)
67 integer :: n
68
69 namelist /fld_list/ fld_kappa0, fld_radforce_split, &
74
75 do n = 1, size(files)
76 open(unitpar, file=trim(files(n)), status="old")
77 read(unitpar, fld_list, end=111)
78 111 close(unitpar)
79 end do
80 end subroutine fld_params_read
81
82 !> Initialising FLD-module
83 !> Read opacities
84 !> Initialise Multigrid and adimensionalise kappa
85 subroutine fld_init(r_gamma)
88 use mod_physics
91
92 double precision, intent(in) :: r_gamma
93
95 fld_debug=.false.
96 fld_no_mg=.false.
97 ! initialize constant opacity with free electron Thomson scattering value
100 ! sanity checks on input
101 if(fld_kappa0<smalldouble)then
102 if(mype==0) print *,'fld_kappa0=',fld_kappa0
103 call mpistop("please set the constant opacity to a reasonable value")
104 endif
105 if(fld_bisect_tol<smalldouble)then
106 if(mype==0) print *,'fld_bisect_tol=',fld_bisect_tol
107 call mpistop("convergence tolerance for root solver too strict")
108 endif
109 if(fld_tiring_explicit)then
110 if(mype==0) print *,'Will do photon tiring explicit!!!'
111 endif
112 if(.not.fld_no_mg)then
113 if(fld_diff_tol<smalldouble)then
114 if(mype==0) print *,'fld_diff_tol=',fld_diff_tol
115 call mpistop("convergence tolerance for MG solver too strict")
116 endif
117 select case(nth_for_diff_mg)
118 case(1)
119 ! no need for stencil extension
120 case(2)
121 ! need for stencil extension
123 case default
124 call mpistop("nth_for_diff_mg must be 1 or 2")
125 end select
129 ! store the diffusion coefficient as extra variable (needed in mg vhelmholtz)
130 i_diff_mg = var_set_extravar("D", "D")
131 use_multigrid = .true.
132 ! use multigrid to solve a helmholtz equation with variable coefficient
133 ! this is stored also as extra variable in the mg solver as mg_iveps
134 mg%n_extra_vars = 1
135 mg%operator_type = mg_vhelmholtz
136 ! choice of smoother: can be mg_smoother_gs or gsrb (latter recommended)
137 mg%smoother_type = mg_smoother_gsrb
138 endif
139 !> set gamma
140 fld_gamma = r_gamma
141 !> Read in opacity table if necesary
142 if(trim(fld_opacity_law) .eq. 'opal') then
143 if(si_unit)call mpistop("adjust opal module with SI-cgs conversions for SI - or use cgs!")
145 endif
146 end subroutine fld_init
147
148 !> Set the boundaries for the diffusion of E
153
154 integer :: iB
155
156 ! Set boundary conditions for the multigrid solver
157 do ib = 1, 2*ndim
158 select case (typeboundary(iw_r_e, ib))
159 case (bc_symm)
160 ! d/dx u = 0
161 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
162 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
163 case (bc_asymm)
164 ! u = 0
165 mg%bc(ib, mg_iphi)%bc_type = mg_bc_dirichlet
166 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
167 case (bc_cont)
168 ! d/dx u = 0
169 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
170 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
171 case (bc_periodic)
172 ! Nothing to do here, this is picked up through periodB variable
173 case (bc_special)
174 if(mype==0)then
175 print *,'Special boundary for Erad needs specific user-set MG BC treatment'
176 print *,' and this could be through usr_special_mg_bc call'
177 endif
178 if (associated(usr_special_mg_bc)) then
179 call usr_special_mg_bc(ib)
180 endif
181 case default
182 call mpistop("divE_multigrid warning: unknown b.c. ")
183 end select
184 ! Neumann on diffusion coefficient is needed on all lower grid levels
185 ! d/dx u = 0
186 mg%bc(ib, mg_iveps)%bc_type = mg_bc_neumann
187 mg%bc(ib, mg_iveps)%bc_value = 0.0_dp
188 end do
189
190 end subroutine fld_set_mg_bounds
191
192
193 !> w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
194 !> This subroutine handles the radiation force and its work added explicitly
195 !> and the energy interaction term combined with photon tiring using an implicit update
196 subroutine add_fld_rad_force(qdt,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
197 use mod_constants
199 use mod_geometry
201 integer, intent(in) :: ixi^l, ixo^l
202 double precision, intent(in) :: qdt, x(ixi^s,1:ndim)
203 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw)
204 double precision, intent(inout) :: w(ixi^s,1:nw)
205 logical, intent(in) :: qsourcesplit
206 logical, intent(inout) :: active
207
208 integer :: idir,jdir,nth_for_fld,ix^d
209 double precision, dimension(ixI^S) :: a1,a2,a3,c0,c1,kappa
210 double precision, dimension(ixI^S) :: e_gas,e_rad,tmp
211 double precision, dimension(ixI^S,1:ndim,1:ndim) :: div_v,edd
212
213 !> Calculate and add sourceterms
214 if(qsourcesplit .eqv. fld_radforce_split) then
215 active = .true.
216 nth_for_fld=2
217 ! store here lambda in a1 and fld_R in a2
218 call fld_get_eddington(wctprim,x,ixi^l,ixo^l,edd,a1,a2,nth_for_fld)
219
220 !> Photon tiring : calculate tensor grad v (named div_v here)
221 ! NOTE: This is ok for uniform Cartesian only!!!!!
222 ! TODO: introduce gradient of vector in geometry module and call that one
223 do idir = 1,ndim
224 do jdir = 1,ndim
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)
227 enddo
228 enddo
229 ! perform contraction fe : grad(v) with fe eddington tensor
230 {^ifoned
231 a3(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1)
232 }
233 {^iftwod
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)
238 }
239 {^ifthreed
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)
249 }
250
251 do idir = 1,ndim
252 call gradient(wctprim(ixi^s,iw_r_e),ixi^l,ixo^l,idir,tmp,nth_for_fld)
253 ! Radiation force = kappa*rho/c *Flux = lambda gradE
254 ! recycle grad E to store -lambda (grad E)_i
255 tmp(ixo^s) = -a1(ixo^s)*tmp(ixo^s)
256 !> Momentum equation source term
257 w(ixo^s,iw_mom(idir)) = w(ixo^s,iw_mom(idir))+ qdt*tmp(ixo^s)
258 !> Energy equation source term
259 w(ixo^s,iw_e) = w(ixo^s,iw_e) + qdt*wctprim(ixo^s,iw_mom(idir))*tmp(ixo^s)
260 enddo
261 !> photon tiring when handled explicitly
262 if(fld_tiring_explicit)then
263 w(ixo^s,iw_r_e) = w(ixo^s,iw_r_e) - qdt*wctprim(ixo^s,iw_r_e)*a3(ixo^s)
264 a3(ixo^s)=0.0d0
265 endif
266
267 call get_and_check_egas_erad_from_conserved(w,ixi^l,ixo^l,e_gas,e_rad)
268 ! BEGIN radiative exchange part with or without photon tiring included
269 call fld_get_opacity_prim(wctprim,x,ixi^l,ixo^l,kappa)
270 !> Coefficients for the polynomial in Moens et al. 2022, eq 37. but with photon tiring (a3)
271 ! NOTE: the next two lines are to be updated when generic EOS in place
272 call phys_get_rfactor(wct,x,ixi^l,ixo^l,tmp)
273 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)
274 a2(ixo^s) = c_norm*kappa(ixo^s)*wct(ixo^s,iw_rho)*qdt
275 a3(ixo^s) = a3(ixo^s)*qdt
276
277 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))
278 c1(ixo^s) = (one+a2(ixo^s)+a3(ixo^s))/a1(ixo^s)/(one+a3(ixo^s))
279
280 !> Loop over every cell for rootfinding method
281 {do ix^d = ixomin^d,ixomax^d\}
282 select case(fld_interaction_method)
283 case('Bisect')
284 call bisection_method(e_gas(ix^d),c0(ix^d),c1(ix^d))
285 case('Newton')
286 call newton_method(e_gas(ix^d),c0(ix^d),c1(ix^d))
287 case('Halley')
288 call halley_method(e_gas(ix^d),c0(ix^d),c1(ix^d))
289 case default
290 call mpistop('root-method not known')
291 end select
292 {enddo\}
293
294 e_rad(ixo^s) = (a1(ixo^s)*e_gas(ixo^s)**4.d0+e_rad(ixo^s))/(one+a2(ixo^s)+a3(ixo^s))
295
296 if(check_small_values.and..not.fix_small_values)then
297 {do ix^db= ixomin^db,ixomax^db\}
298 if(e_gas(ix^d)<small_e.or.e_rad(ix^d)<small_r_e) then
299 write(*,*) "Error in FLD add_fld_rad_force: small value"
300 write(*,*) "of internal or radiation energy density after exchange"
301 write(*,*) "Iteration: ", it, " Time: ", global_time
302 write(*,*) "Location: ", x(ix^d,:)
303 write(*,*) "Cell number: ", ix^d
304 write(*,*) "internal energy density is=",e_gas(ix^d)," versus small_e=",small_e
305 write(*,*) "radiation energy density is=",e_rad(ix^d)," versus small_r_e=",small_r_e
306 call mpistop("FLD error:May need to turn on fixes")
307 end if
308 {end do\}
309 endif
310
311 if(fix_small_values)then
312 {do ix^d = ixomin^d,ixomax^d\ }
313 e_gas(ix^d) = max(e_gas(ix^d),small_e)
314 e_rad(ix^d) = max(e_rad(ix^d),small_r_e)
315 {enddo\}
316 endif
317
318 ! END radiative exchange part with or without photon tiring included
319
320 !> Update gas-energy in w, internal + kinetic
321 w(ixo^s,iw_e) = e_gas(ixo^s)
322 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)
323 if(allocated(iw_mag)) then
324 w(ixo^s,iw_e) = w(ixo^s,iw_e)+half*sum(w(ixo^s,iw_mag(:))**2,dim=ndim+1)
325 endif
326 !> Update rad-energy in w
327 w(ixo^s,iw_r_e) = e_rad(ixo^s)
328 end if
329 end subroutine add_fld_rad_force
330
331 !> get dt limit for radiation force and FLD explicit source additions
332 !> NOTE: w is primitive on entry
333 subroutine fld_radforce_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
336 use mod_geometry
338
339 integer, intent(in) :: ixi^l, ixo^l
340 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim), w(ixi^s,1:nw)
341 double precision, intent(inout) :: dtnew
342
343 integer :: idim,idims,nth_for_fld
344 double precision :: dxinv(1:ndim), max_grav
345 double precision :: lambda(ixi^s),fld_r(ixi^s)
346 double precision :: tmp(ixi^s),boostfactor
347 double precision :: max_diff,max_diff_coef,max_diff_dim,dtdifflimit
348 double precision :: cmax(ixi^s),cmaxtot(ixi^s),courantmaxtots
349
350 if(fld_debug)print *,'DT limit on entry to radforce_get_dt=',dtnew
351 nth_for_fld=2
352 call fld_get_fluxlimiter_prim(w,x,ixi^l,ixo^l,lambda,fld_r,nth_for_fld)
353 if(slab_uniform) then
354 ^d&dxinv(^d)=one/dx^d;
355 do idim = 1, ndim
356 call gradient(w(ixi^s,iw_r_e),ixi^l,ixo^l,idim,tmp,nth_for_fld)
357 max_grav = maxval(dabs(-lambda(ixo^s)*tmp(ixo^s)/w(ixo^s,iw_rho)))
358 max_grav = max(max_grav, epsilon(1.0d0))
359 dtnew = min(dtnew, 1.0d0 / dsqrt(max_grav * dxinv(idim)))
360 end do
361 else
362 do idim = 1, ndim
363 call gradient(w(ixi^s,iw_r_e),ixi^l,ixo^l,idim,tmp,nth_for_fld)
364 max_grav = maxval(dabs(-lambda(ixo^s)*tmp(ixo^s)/w(ixo^s,iw_rho))/block%ds(ixo^s,idim))
365 max_grav = max(max_grav, epsilon(1.0d0))
366 dtnew = min(dtnew, 1.0d0 / dsqrt(max_grav))
367 end do
368 endif
369 if(fld_debug)print *,'DT limit after RADFORCE eff grav=',dtnew
370
371 if(fld_slowsteps) then
372 boostfactor=1.0d0
373 call fld_get_opacity_prim(w, x, ixi^l, ixo^l, tmp)
374 max_diff=0.0d0
375 if(slab_uniform) then
376 max_diff_coef = maxval(c_norm*lambda(ixo^s)/(w(ixo^s,iw_rho)*tmp(ixo^s)))
377 ^d&dxinv(^d)=one/dx^d;
378 do idim = 1, ndim
379 max_diff_dim = max(2.0d0*ndim*max_diff_coef*dxinv(idim)**2, epsilon(1.0d0))
380 max_diff = max(max_diff,max_diff_dim)
381 end do
382 else
383 do idim = 1, ndim
384 max_diff_coef = maxval(c_norm*lambda(ixo^s)/(w(ixo^s,iw_rho)*tmp(ixo^s))/block%ds(ixo^s,idim)**2)
385 max_diff_dim = max(2.0d0*ndim*max_diff_coef, epsilon(1.0d0))
386 max_diff = max(max_diff,max_diff_dim)
387 end do
388 endif
389 dtdifflimit=1.0d0/max_diff
390 if(fld_debug) print *,'DT limit from diffusion is actually=',dtdifflimit
391 if(slowsteps>it-it_init+1) then
392 boostfactor=1.0d0+fld_boost_dt*(dble(it-it_init+1)/dble(slowsteps))
393 else
394 boostfactor=1.0d0+fld_boost_dt
395 endif
396 dtdifflimit=dtdifflimit*boostfactor
397 if(fld_debug) print *,'DT limit from diffusion is boosted to=',dtdifflimit
398 dtnew=min(dtnew,dtdifflimit)
399 if(fld_debug) print *,'DT limit after diffusion is =',dtnew
400 endif
401
402 ! here we interface back to fld_get_radpress
403 call phys_get_csrad2(w,x,ixi^l,ixo^l,tmp)
404 if(slab_uniform) then
405 ^d&dxinv(^d)=one/dx^d;
406 do idims=1,ndim
407 cmax(ixo^s)=dabs(w(ixo^s,iw_mom(idims)))+dsqrt(tmp(ixo^s))
408 if(idims==1) then
409 cmaxtot(ixo^s)=cmax(ixo^s)*dxinv(idims)
410 else
411 cmaxtot(ixo^s)=cmaxtot(ixo^s)+cmax(ixo^s)*dxinv(idims)
412 end if
413 end do
414 else
415 do idims=1,ndim
416 cmax(ixo^s)=dabs(w(ixo^s,iw_mom(idims)))+dsqrt(tmp(ixo^s))
417 if(idims==1) then
418 cmaxtot(ixo^s)=cmax(ixo^s)/block%ds(ixo^s,idims)
419 else
420 cmaxtot(ixo^s)=cmaxtot(ixo^s)+cmax(ixo^s)/block%ds(ixo^s,idims)
421 end if
422 end do
423 end if
424 ! courantmaxtots='max(summed c/dx)'
425 courantmaxtots=maxval(cmaxtot(ixo^s))
426 if(fld_debug)print *,'DT limit RADFORCE CSRAD=',courantpar/courantmaxtots
427 if(courantmaxtots>smalldouble) dtnew=min(dtnew,courantpar/courantmaxtots)
428 if(fld_debug)print *,'DT limit FINALLY ENFORCED IS NOW=',dtnew
429
430 end subroutine fld_radforce_get_dt
431
432 !> Sets the opacity in the w-array
433 !> by calling mod_opal_opacity
434 !> NOTE: assumes primitives in w
435 !> NOTE: assuming opacity is local, not ok with cak line force
436 subroutine fld_get_opacity_prim(w, x, ixI^L, ixO^L, fld_kappa)
438 use mod_physics, only: phys_get_tgas
441 integer, intent(in) :: ixi^l, ixo^l
442 double precision, intent(in) :: w(ixi^s, 1:nw)
443 double precision, intent(in) :: x(ixi^s, 1:ndim)
444 double precision, intent(out) :: fld_kappa(ixi^s)
445
446 integer :: ix^d
447 double precision :: rho0,temp0,kapp0
448 double precision :: temp(ixi^s)
449
450 select case (trim(fld_opacity_law))
451 case('const_norm')
452 fld_kappa(ixo^s) = fld_kappa0
453 case('const')
454 fld_kappa(ixo^s) = fld_kappa0/unit_opacity
455 case('opal')
456 call phys_get_tgas(w,x,ixi^l,ixo^l,temp)
457 {do ix^d=ixomin^d,ixomax^d\ }
458 rho0 = w(ix^d,iw_rho)*unit_density
459 temp0 = temp(ix^d)*unit_temperature
460 call set_opal_opacity(rho0,temp0,kapp0)
461 fld_kappa(ix^d) = kapp0/unit_opacity
462 {enddo\ }
463 case('special')
464 if (.not. associated(usr_special_opacity)) then
465 call mpistop("special opacity not defined")
466 endif
467 call usr_special_opacity(ixi^l, ixo^l, w, x, fld_kappa)
468 case default
469 call mpistop("Doesn't know opacity law")
470 end select
471 end subroutine fld_get_opacity_prim
472
473 !> Returns Radiation Pressure as tensor
474 !> NOTE: w is primitive on entry
475 subroutine fld_get_radpress(w, x, ixI^L, ixO^L, rad_pressure)
477 integer, intent(in) :: ixi^l, ixo^l
478 double precision, intent(in) :: w(ixi^s, 1:nw)
479 double precision, intent(in) :: x(ixi^s, 1:ndim)
480 double precision, intent(out):: rad_pressure(ixi^s,1:ndim,1:ndim)
481 integer :: i,j,nth
482 double precision :: eddington_tensor(ixi^s,1:ndim,1:ndim)
483 double precision :: lambda(ixi^s),fld_r(ixi^s)
484
485 ! always use 4th order CD here
486 nth=2
487 call fld_get_eddington(w, x, ixi^l, ixo^l, eddington_tensor, lambda, fld_r, nth)
488 if(fld_debug.and..false.)then
489 print *,'In get_radPress with nth=',nth,' on ixO=',ixo^l
490 print *,'Max and Min value of fe'
491 print *,maxval(eddington_tensor(ixo^s,1:ndim,1:ndim))
492 print *,minval(eddington_tensor(ixo^s,1:ndim,1:ndim))
493 print *,'Max and Min value of Erad'
494 print *,maxval(w(ixo^s,iw_r_e))
495 print *,minval(w(ixo^s,iw_r_e))
496 print *,'End get_radPress'
497 endif
498 do i=1,ndim
499 do j=1,ndim
500 rad_pressure(ixo^s,i,j) = eddington_tensor(ixo^s,i,j)*w(ixo^s,iw_r_e)
501 enddo
502 enddo
503 end subroutine fld_get_radpress
504
505 !> This subroutine calculates flux limiter lambda according to fld_fluxlimiter
506 !> It also calculates fld_R which is ratio of radiation scaleheight and mean free path
507 !> NOTE: nth and ixI and ixO not free to choose here: TODO
508 subroutine fld_get_fluxlimiter(w,x,ixI^L,ixO^L,fld_lambda,fld_R,nth)
510 use mod_geometry
513 integer, intent(in) :: ixi^l,ixo^l,nth
514 double precision, intent(in) :: w(ixi^s,1:nw)
515 double precision, intent(in) :: x(ixi^s,1:ndim)
516 double precision, intent(out) :: fld_r(ixi^s),fld_lambda(ixi^s)
517
518 double precision :: wprim(ixi^s,1:nw)
519
520 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
521 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
522 call fld_get_fluxlimiter_prim(wprim,x,ixi^l,ixo^l,fld_lambda,fld_r,nth)
523
524 end subroutine fld_get_fluxlimiter
525
526 !> This subroutine calculates flux limiter lambda according to fld_fluxlimiter
527 !> It also calculates fld_R which is ratio of radiation scaleheight and mean free path
528 !> NOTE: this one operates on primitives
529 !> NOTE: nth and ixI and ixO not free to choose
530 subroutine fld_get_fluxlimiter_prim(w,x,ixI^L,ixO^L,fld_lambda,fld_R,nth)
532 use mod_geometry
534 integer, intent(in) :: ixi^l,ixo^l,nth
535 double precision, intent(in) :: w(ixi^s,1:nw)
536 double precision, intent(in) :: x(ixi^s,1:ndim)
537 double precision, intent(out) :: fld_r(ixi^s),fld_lambda(ixi^s)
538
539 integer :: idir, ix^d
540 double precision :: kappa(ixi^s),normgrad2(ixi^s)
541 double precision :: grad_r_e(ixi^s)
542
543 select case(fld_fluxlimiter)
544 case('Diffusion')
545 ! optically thick limit
546 fld_lambda(ixo^s) = 1.d0/3.d0
547 fld_r(ixo^s) = zero
548 case('FreeStream')
549 ! optically thin limit
550 normgrad2(ixo^s) = zero
551 do idir=1,ndim
552 call gradient(w(ixi^s,iw_r_e),ixi^l,ixo^l,idir,grad_r_e,nth)
553 normgrad2(ixo^s) = normgrad2(ixo^s)+grad_r_e(ixo^s)**2
554 end do
555 call fld_get_opacity_prim(w,x,ixi^l,ixo^l,kappa)
556 ! Calculate R everywhere
557 ! |grad E|/(rho kappa E)
558 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
559 where(normgrad2(ixo^s)<smalldouble**2)
560 ! treat uniform case as diffusion limit
561 fld_r(ixo^s)=zero
562 fld_lambda(ixo^s) = 1.0d0/3.0d0
563 elsewhere
564 fld_lambda(ixo^s) = one/fld_r(ixo^s)
565 endwhere
566 case('Pomraning')
567 ! Calculate R everywhere
568 ! |grad E|/(rho kappa E)
569 normgrad2(ixo^s) = zero
570 do idir = 1,ndim
571 call gradient(w(ixi^s,iw_r_e),ixi^l,ixo^l,idir,grad_r_e,nth)
572 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_e(ixo^s)**2
573 end do
574 call fld_get_opacity_prim(w,x,ixi^l,ixo^l,kappa)
575 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
576 ! Calculate the flux limiter, lambda
577 ! Levermore and Pomraning: lambda = (2 + R)/(6 + 3R + R^2)
578 fld_lambda(ixo^s) = (2.d0+fld_r(ixo^s))/(6.d0+3*fld_r(ixo^s)+fld_r(ixo^s)**2)
579 case('Minerbo')
580 ! Calculate R everywhere
581 ! |grad E|/(rho kappa E)
582 normgrad2(ixo^s) = zero
583 do idir = 1,ndim
584 call gradient(w(ixi^s,iw_r_e),ixi^l,ixo^l,idir,grad_r_e,nth)
585 normgrad2(ixo^s) = normgrad2(ixo^s) + grad_r_e(ixo^s)**2
586 end do
587 call fld_get_opacity_prim(w, x, ixi^l, ixo^l, kappa)
588 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
589 ! Calculate the flux limiter, lambda
590 ! Minerbo:
591 {do ix^d = ixomin^d,ixomax^d\ }
592 if(fld_r(ix^d) .lt. 3.d0/2.d0) then
593 fld_lambda(ix^d) = 2.d0/(3.d0+dsqrt(9.d0+12.d0*fld_r(ix^d)**2))
594 else
595 fld_lambda(ix^d) = 1.d0/(1.d0+fld_r(ix^d)+dsqrt(1.d0+2.d0*fld_r(ix^d)))
596 endif
597 {enddo\}
598 case('special')
599 if (.not. associated(usr_special_fluxlimiter)) then
600 call mpistop("special fluxlimiter not defined")
601 endif
602 call usr_special_fluxlimiter(ixi^l,ixo^l,w,x,fld_lambda,fld_r)
603 case default
604 call mpistop('Fluxlimiter unknown')
605 end select
606
607 end subroutine fld_get_fluxlimiter_prim
608
609 !> Calculate Radiation Flux
610 !> NOTE: only for diagnostics purposes (w conservative on entry)
611 !> This returns cell centered values for radiation flux
612 subroutine fld_get_radflux(w, x, ixI^L, ixO^L, rad_flux)
614 use mod_geometry
616 integer, intent(in) :: ixi^l, ixo^l
617 double precision, intent(in) :: w(ixi^s, 1:nw)
618 double precision, intent(in) :: x(ixi^s, 1:ndim)
619 double precision, intent(out) :: rad_flux(ixi^s, 1:ndim)
620
621 integer :: idir,nth_for_fld
622 double precision :: wprim(ixi^s,1:nw)
623 double precision :: grad_r_e(ixi^s)
624 double precision :: kappa(ixi^s), lambda(ixi^s), fld_r(ixi^s)
625
626 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
627 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
628
629 call fld_get_opacity_prim(wprim, x, ixi^l, ixo^l, kappa)
630 ! always use 4th order CD here
631 nth_for_fld=2
632 call fld_get_fluxlimiter_prim(wprim, x, ixi^l, ixo^l, lambda, fld_r, nth_for_fld)
633 !> Calculate the Flux using the fld closure relation
634 !> F = -c*lambda/(kappa*rho) *grad E
635 do idir = 1,ndim
636 call gradient(wprim(ixi^s,iw_r_e),ixi^l,ixo^l,idir,grad_r_e,nth_for_fld)
637 rad_flux(ixo^s,idir)=-(c_norm*lambda(ixo^s)/(kappa(ixo^s)*wprim(ixo^s,iw_rho)))*grad_r_e(ixo^s)
638 end do
639 end subroutine fld_get_radflux
640
641 !> Calculate Eddington-tensor (where w is primitive)
642 !> also feeds back the flux limiter lambda and R
643 subroutine fld_get_eddington(w, x, ixI^L, ixO^L, eddington_tensor, lambda, fld_R, nth)
645 use mod_geometry
646 integer, intent(in) :: ixI^L, ixO^L, nth
647 double precision, intent(in) :: w(ixI^S, 1:nw)
648 double precision, intent(in) :: x(ixI^S, 1:ndim)
649 double precision, intent(out) :: eddington_tensor(ixI^S,1:ndim,1:ndim)
650 double precision, intent(out) :: lambda(ixI^S),fld_R(ixI^S)
651
652 integer :: idir,jdir
653 double precision :: normgrad2(ixI^S)
654 double precision :: tmp(ixI^S),grad_r_e(ixI^S,1:ndim)
655 double precision :: nn_regularized(ixI^S,1:ndim,1:ndim)
656
657 normgrad2(ixo^s) = zero
658 do idir = 1,ndim
659 call gradient(w(ixi^s, iw_r_e),ixi^l,ixo^l,idir,tmp,nth)
660 grad_r_e(ixo^s,idir)=tmp(ixo^s)
661 normgrad2(ixo^s)=normgrad2(ixo^s)+tmp(ixo^s)**2
662 end do
663 do idir = 1,ndim
664 do jdir = 1,ndim
665 if(idir==jdir)then
666 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)
667 else
668 nn_regularized(ixo^s,idir,jdir)=(grad_r_e(ixo^s,idir)*grad_r_e(ixo^s,jdir))/(normgrad2(ixo^s)+smalldouble**2)
669 endif
670 enddo
671 enddo
672 ! get lambda and R
673 call fld_get_fluxlimiter_prim(w,x,ixi^l,ixo^l,lambda,fld_r,nth)
674 ! store f_e= lambda + lambda^2 R^2
675 tmp(ixo^s) = lambda(ixo^s)+(lambda(ixo^s)*fld_r(ixo^s))**2
676 do idir = 1,ndim
677 ! first compute the isotropic (diagonal) part
678 eddington_tensor(ixo^s,idir,idir) = half*(one-tmp(ixo^s))
679 enddo
680 do idir = 1,ndim
681 do jdir = 1,ndim
682 ! initialize off-diagonal part here
683 if(idir .ne. jdir) eddington_tensor(ixo^s,idir,jdir) = zero
684 ! add part depending on unit vectors along gradient E
685 eddington_tensor(ixo^s,idir,jdir) = eddington_tensor(ixo^s,idir,jdir)+&
686 half*(3.d0*tmp(ixo^s)-one)*nn_regularized(ixo^s,idir,jdir)
687 enddo
688 enddo
689 end subroutine fld_get_eddington
690
691 subroutine get_and_check_egas_erad_from_conserved(w,ixI^L,ixO^L,e_gas,E_rad)
693 integer, intent(in) :: ixI^L, ixO^L
694 double precision, intent(in) :: w(ixI^S, 1:nw)
695 double precision, intent(out) :: e_gas(ixI^S),E_rad(ixI^S)
696 integer :: ix^D
697
698 !> e_gas is the INTERNAL ENERGY without KINETIC ENERGY
699 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)
700 if(allocated(iw_mag)) then
701 e_gas(ixo^s) = e_gas(ixo^s)-half*sum(w(ixo^s,iw_mag(:))**2,dim=ndim+1)
702 endif
703 e_rad(ixo^s) = w(ixo^s,iw_r_e)
704
705 if(check_small_values.and..not.fix_small_values)then
706 {do ix^db= ixomin^db,ixomax^db\}
707 if(e_gas(ix^d)<small_e.or.e_rad(ix^d)<small_r_e) then
708 write(*,*) "Error in FLD get_egas_Erad: small value"
709 write(*,*) "Iteration: ", it, " Time: ", global_time
710 write(*,*) "Cell number: ", ix^d
711 write(*,*) "internal energy density is=",e_gas(ix^d)," versus small_e=",small_e
712 write(*,*) "radiation energy density is=",e_rad(ix^d)," versus small_r_e=",small_r_e
713 call mpistop("FLD error:May need to turn on fixes")
714 end if
715 {end do\}
716 endif
717
718 if(fix_small_values)then
719 {do ix^d = ixomin^d,ixomax^d\ }
720 e_gas(ix^d) = max(e_gas(ix^d),small_e)
721 e_rad(ix^d) = max(e_rad(ix^d),small_r_e)
722 {enddo\}
723 endif
724
726
727 !> Calling all subroutines to perform the multigrid method
728 !> Communicates rad_e and diff_coeff to multigrid library
729 !> Advance psa=psb+dtfactor*qdt*F_im(psa)
730 subroutine fld_implicit_update(dtfactor,qdt,qtC,psa,psb)
732 use mod_forest
735
736 type(state), target :: psa(max_blocks)
737 type(state), target :: psb(max_blocks)
738 double precision, intent(in) :: qdt
739 double precision, intent(in) :: qtC
740 double precision, intent(in) :: dtfactor
741
742 integer, parameter :: max_its = 100
743 integer :: n,ixO^L,ix^D
744 double precision :: res, max_residual, mg_lambda, fac
745 double precision :: wmax(nw),wmin(nw)
746 double precision :: wmaxb(nw),wminb(nw)
747 integer :: iigrid, igrid
748
749 ! we need first to compute the (variable) diffusion coefficient on entire grid
750 ! this must be done in mesh+1 ghostcell layer
751 call update_diffcoeff(psa)
752
753 ! now we multiply the diffusion coefficient with dtfactor*dt on entire mesh+1 domain
754 ixo^l=ixm^ll^ladd1;
755 do iigrid=1,igridstail; igrid=igrids(iigrid);
756 {do ix^d = ixomin^d,ixomax^d\ }
757 psa(igrid)%w(ix^d,i_diff_mg)= psa(igrid)%w(ix^d,i_diff_mg)*dtfactor*qdt
758 {enddo\}
759 end do
760
761 fac = 1.0d0
762 max_residual = fld_diff_tol
763
764 if(fld_debug)then
765 call get_global_maxima(wmax,psa)
766 call get_global_minima(wmin,psa)
767 call get_global_maxima(wmaxb,psb)
768 call get_global_minima(wminb,psb)
769 if(mype==0)then
770 ! the MG needs to be scaled such that everything is order unity
771 print *,'Currently at time=',global_time,' time step=',qdt,' dtfactor=',dtfactor
772 print *,'at start of MG solver, we have fld_diff_tol =',fld_diff_tol
773 print *,'at start of MG solver, we have LHS E_rad range as :',wmax(iw_r_e),wmin(iw_r_e)
774 print *,'at start of MG solver, we have Diff coeff range as:',wmax(i_diff_mg),wmin(i_diff_mg)
775 print *,'at start of MG solver, we have density range as :',wmax(iw_rho),wmin(iw_rho)
776 print *,'at start of MG solver, we have RHS E_rad range as :',wmaxb(iw_r_e),wminb(iw_r_e)
777 print *,'at start of MG solver, we have qdt as',qdt,' and max_residual=',max_residual
778 print *,'at start of MG solver, ratio coeffs on level 1 =',wmax(i_diff_mg)/(dx(1,1)**2)
779 print *,'at start of MG solver, ratio coeffs on level max=',wmax(i_diff_mg)/(dx(1,refine_max_level)**2)
780 endif
781 endif
782
783 call mg_set_methods(mg)
784 if(.not. mg%is_allocated) call mpistop("multigrid tree not allocated yet")
785
786 ! Here we handle the global helmholtz problem with variable coefficient
787 ! The equation we solve is div([D]^n nabla Erad^(n+1)) -(1/dt)Erad^(n+1)=-(1/dt)Erad^n
788 ! we reformulate to div([D x dt]^n nabla Erad^(n+1)) -Erad^(n+1)=-Erad^n
789 ! Helmholtz equation is div(eps nabla phi) -mg_lambda phi = f
790 ! hence eps is our variable coefficient dt x D=fld_lambda*c/(kappa*rho)
791 ! hence phi is Erad and mg_lambda is unity
792 call vhelmholtz_set_lambda(fac)
793 ! copy in the (variable) diffusion coefficient in mg_iveps
794 ! NOTE: this copies also the ghostcell values for this coefficient to mg
795 call mg_copy_to_tree(i_diff_mg, mg_iveps, factor=fac, state_from=psa)
796 ! copy in the Erad variable in mg_iphi (the one we solve for)
797 call mg_copy_to_tree(iw_r_e, mg_iphi, factor=fac, state_from=psa)
798 ! copy in RHS f factor as Erad with factor -1
799 call mg_copy_to_tree(iw_r_e, mg_irhs, factor=-fac, state_from=psb)
800 ! becuase the variable coefficient is needed on lower grid levels in mg
801 ! we need to restrict and adopt BCs for this variable: Neumann is set
802 call mg_restrict(mg, mg_iveps)
803 call mg_fill_ghost_cells(mg, mg_iveps)
804 ! Now try solving with MG
805 call mg_fas_fmg(mg, .true., max_res=res)
806 if(fld_debug.and.mype==0)print *,'MG residual obtained with FMG is =',res
807 do n = 1, max_its
808 if(res < max_residual) exit
809 call mg_fas_vcycle(mg, max_res=res)
810 if(fld_debug.and.mype==0)print *,'MG residual obtained with Vcycle =',res,' at step =',n
811 end do
812 if(fld_debug.and.mype==0)print *,'FINAL MG residual obtained is =',res
813 if(res .ge. max_residual) then
814 if (mype == 0) then
815 write(*,*) it, ' residual from MG ', res
816 write(*,*) it, ' max_residual in MG ', max_residual
817 write(*,*) it, ' dtfactor*qdt in MG ', qdt*dtfactor
818 print *,'Currently at time=',global_time,' time step=',qdt,' dtfactor=',dtfactor
819 endif
821 call mpistop("no convergence in MG")
822 else
823 if(mype==0) write(*,*) 'WARNING for it=',it,' NO CONVERGENCE IN MG but still carry on'
824 endif
825 end if
826 ! copy back the Erad variable in iw_r_e
827 call mg_copy_from_tree_gc(mg_iphi, iw_r_e, state_to=psa)
828
829 if(check_small_values.and..not.fix_small_values)then
830 ixo^l=ixm^ll;
831 do iigrid=1,igridstail; igrid=igrids(iigrid);
832 {do ix^db= ixomin^db,ixomax^db\}
833 if(psa(igrid)%w(ix^d,iw_r_e)<small_r_e) then
834 write(*,*) "Error in FLD fld_implicit_update: small value"
835 write(*,*) "of radiation energy density after MG"
836 write(*,*) "Iteration: ", it, " Time: ", global_time
837 write(*,*) "Location: ", psa(igrid)%x(ix^d,:)," on grid",igrid
838 write(*,*) "Cell number: ", ix^d
839 write(*,*) "radiation energy density is=",psa(igrid)%w(ix^d,iw_r_e)," versus small_r_e=",small_r_e
840 call mpistop("FLD error:May need to turn on fixes")
841 end if
842 {end do\}
843 end do
844 endif
845
846 if(fix_small_values)then
847 ixo^l=ixm^ll;
848 do iigrid=1,igridstail; igrid=igrids(iigrid);
849 {do ix^d = ixomin^d,ixomax^d\ }
850 psa(igrid)%w(ix^d,iw_r_e)= max(psa(igrid)%w(ix^d,iw_r_e),small_r_e)
851 {enddo\}
852 end do
853 endif
854
855 end subroutine fld_implicit_update
856
857 subroutine update_diffcoeff(psa)
859 type(state), target :: psa(max_blocks)
860
861 integer :: iigrid, igrid, ixO^L
862
863 ! we will need diffusion coefficients in i-1 i i+1
864 ixo^l=ixm^ll^ladd1;
865 !$OMP PARALLEL DO PRIVATE(igrid)
866 do iigrid=1,igridstail; igrid=igrids(iigrid);
867 call fld_get_diffcoef_central(psa(igrid)%w, psa(igrid)%x, ixg^ll, ixo^l)
868 end do
869 !$OMP END PARALLEL DO
870
871 end subroutine update_diffcoeff
872
873 !> inplace update of psa==>F_im(psa)
874 subroutine fld_evaluate_implicit(qtC,psa)
876 type(state), target :: psa(max_blocks)
877 double precision, intent(in) :: qtC
878 integer :: iigrid, igrid
879 integer :: ixO^L
880
881 call update_diffcoeff(psa)
882
883 ixo^l=ixm^ll;
884 !$OMP PARALLEL DO PRIVATE(igrid)
885 do iigrid=1,igridstail; igrid=igrids(iigrid);
886 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
887 call evaluate_diffterm_onegrid(ixg^ll,ixo^l,psa(igrid)%w,psa(igrid)%x)
888 end do
889 !$OMP END PARALLEL DO
890
891 end subroutine fld_evaluate_implicit
892
893 !> inplace update of psa==>F_im(psa)
894 subroutine evaluate_diffterm_onegrid(ixI^L,ixO^L,w,x)
896 integer, intent(in) :: ixI^L, ixO^L
897 double precision, intent(inout) :: w(ixI^S, 1:nw)
898 double precision, intent(in) :: x(ixI^S, 1:ndim)
899
900 double precision :: divF(ixI^S)
901 integer :: idir, jxO^L, hxO^L
902
903 if(.not.slab)call mpistop("laplacian coded up for uniform cartesian grid")
904
905 ! Here we use diffusion coefficients in positions i-1 i i+1 and exploit harmonic means i.e. 2ab/(a+b)
906 ! since this is how the multigrid library handles the div(eps nabla phi) term in Helmholtz equation
907 ! div(eps nabla phi) - mg_lambda phi = f
908 divf(ixo^s) = 0.d0
909 do idir = 1,ndim
910 hxo^l=ixo^l-kr(idir,^d);
911 jxo^l=ixo^l+kr(idir,^d);
912 divf(ixo^s) = divf(ixo^s) + &
913 (w(jxo^s,iw_r_e)*two*w(ixo^s,i_diff_mg)*w(jxo^s,i_diff_mg)/(w(ixo^s,i_diff_mg) + w(jxo^s,i_diff_mg)) &
914 -w(ixo^s,iw_r_e)*(two*w(ixo^s,i_diff_mg)*w(jxo^s,i_diff_mg)/(w(ixo^s,i_diff_mg) + w(jxo^s,i_diff_mg)) &
915 +two*w(ixo^s,i_diff_mg)*w(hxo^s,i_diff_mg)/(w(ixo^s,i_diff_mg) + w(hxo^s,i_diff_mg))) &
916 +w(hxo^s,iw_r_e)*two*w(ixo^s,i_diff_mg)*w(hxo^s,i_diff_mg)/(w(ixo^s,i_diff_mg) + w(hxo^s,i_diff_mg)))/dxlevel(idir)**2
917 ! below uses artihmetic mean, different from mg method
918 !divF(ixO^S) = divF(ixO^S) + &
919 ! (w(jxO^S,iw_r_e)*half*(w(ixO^S,i_diff_mg) + w(jxO^S,i_diff_mg)) &
920 ! -w(ixO^S,iw_r_e)*(half*(w(ixO^S,i_diff_mg) + w(jxO^S,i_diff_mg)) &
921 ! +half*(w(ixO^S,i_diff_mg) + w(hxO^S,i_diff_mg))) &
922 ! +w(hxO^S,iw_r_e)*half*(w(ixO^S,i_diff_mg) + w(hxO^S,i_diff_mg)))/dxlevel(idir)**2
923 enddo
924 ! only the E variable is handled implicitly, all else must be zero here
925 w(ixo^s,1:nw)=zero
926 w(ixo^s,iw_r_e) = divf(ixo^s)
927
928 end subroutine evaluate_diffterm_onegrid
929
930 !> Calculates cell-centered diffusion coefficient to be used in multigrid
931 subroutine fld_get_diffcoef_central(w, x, ixI^L, ixO^L)
933 use mod_geometry
936 integer, intent(in) :: ixi^l, ixo^l
937 double precision, intent(in) :: x(ixi^s,1:ndim)
938 double precision, intent(inout) :: w(ixi^s,1:nw)
939
940 double precision :: wprim(ixi^s,1:nw)
941 double precision :: kappa(ixi^s),lambda(ixi^s),fld_r(ixi^s)
942
943 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
944 ! ensure entries in entire ixI range
945 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
946 call fld_get_opacity_prim(wprim, x, ixi^l, ixo^l, kappa)
947 ! note that we use central difference here (last argument is 1 or 2)
948 call fld_get_fluxlimiter_prim(wprim, x, ixi^l, ixo^l, lambda, fld_r, nth_for_diff_mg)
949 w(ixo^s,i_diff_mg) = c_norm*lambda(ixo^s)/(kappa(ixo^s)*wprim(ixo^s,iw_rho))
950 if(associated(usr_special_diffcoef)) call usr_special_diffcoef(w, wprim, x, ixi^l, ixo^l)
951 if(minval(w(ixo^s,i_diff_mg))<smalldouble) then
952 print *,'min diffcoef=',minval(w(ixo^s,i_diff_mg))
953 call mpistop("too small diffusion coefficient")
954 endif
955 if(maxval(w(ixo^s,i_diff_mg))>bigdouble) call mpistop("too large diffusion coefficient")
956
957 if(fld_debug.and..false.)then
958 print *,'setting diffcoefs with data on',ixi^l
959 print *,'min diffcoef=',minval(w(ixo^s,i_diff_mg))
960 print *,'min lambda kappa rho fld_R'
961 print *,minval(lambda(ixo^s))
962 print *,minval(kappa(ixo^s))
963 print *,minval(wprim(ixo^s,iw_rho))
964 print *,minval(fld_r(ixo^s))
965 print *,'max diffcoef=',maxval(w(ixo^s,i_diff_mg))
966 print *,'max lambda kappa rho fld_R'
967 print *,maxval(lambda(ixo^s))
968 print *,maxval(kappa(ixo^s))
969 print *,maxval(wprim(ixo^s,iw_rho))
970 print *,maxval(fld_r(ixo^s))
971 print *,'done setting diffcoefs in slot',i_diff_mg,' on range',ixo^l
972 endif
973
974 end subroutine fld_get_diffcoef_central
975
976
977 !> Find the root of the 4th degree polynomial using the bisection method
978 subroutine bisection_method(e_gas, c0, c1)
980 double precision, intent(in) :: c0, c1
981 double precision, intent(inout) :: e_gas
982 double precision :: bisect_a, bisect_b, bisect_c
983 integer :: n, max_its
984
985 n = 0
986 max_its = 100
987 bisect_a = zero
988 bisect_b = min(dabs(c0/c1),dabs(c0)**(1.d0/4.d0))+smalldouble
989 do while(dabs(bisect_b-bisect_a) .ge. fld_bisect_tol)
990 bisect_c = (bisect_a + bisect_b)/two
991 n = n +1
992 if(n .gt. max_its) then
993 call mpistop('No convergece in bisection scheme')
994 endif
995 if(polynomial_bisection(bisect_a, c0, c1)*&
996 polynomial_bisection(bisect_b, c0, c1) .lt. zero) then
997 if(polynomial_bisection(bisect_a, c0, c1)*&
998 polynomial_bisection(bisect_c, c0, c1) .lt. zero) then
999 bisect_b = bisect_c
1000 elseif(polynomial_bisection(bisect_b, c0, c1)*&
1001 polynomial_bisection(bisect_c, c0, c1) .lt. zero) then
1002 bisect_a = bisect_c
1003 elseif(polynomial_bisection(bisect_a, c0, c1) .eq. zero) then
1004 bisect_b = bisect_a
1005 bisect_c = bisect_a
1006 goto 2435
1007 elseif(polynomial_bisection(bisect_b, c0, c1) .eq. zero) then
1008 bisect_a = bisect_b
1009 bisect_c = bisect_b
1010 goto 2435
1011 elseif(polynomial_bisection(bisect_c, c0, c1) .eq. zero) then
1012 bisect_a = bisect_c
1013 bisect_b = bisect_c
1014 goto 2435
1015 else
1016 call mpistop("Problem with fld bisection method")
1017 endif
1018 elseif(polynomial_bisection(bisect_a, c0, c1) &
1019 - polynomial_bisection(bisect_b, c0, c1) .lt. fld_bisect_tol*polynomial_bisection(bisect_a, c0, c1)) then
1020 goto 2435
1021 else
1022 bisect_a = e_gas
1023 bisect_b = e_gas
1024 if(fld_debug)print*, "IGNORING GAS-RAD ENERGY EXCHANGE ", c0, c1
1025 if(fld_debug)print*, polynomial_bisection(bisect_a, c0, c1), polynomial_bisection(bisect_b, c0, c1)
1026 call mpistop('issues in bisection scheme')
1027 if(polynomial_bisection(bisect_a, c0, c1) .le. smalldouble) then
1028 bisect_b = bisect_a
1029 elseif(polynomial_bisection(bisect_a, c0, c1) .le. smalldouble) then
1030 bisect_a = bisect_b
1031 endif
1032 goto 2435
1033 endif
1034 enddo
1035 2435 e_gas = (bisect_a + bisect_b)/two
1036 end subroutine bisection_method
1037
1038 !> Find the root of the 4th degree polynomial using the Newton method
1039 subroutine newton_method(e_gas, c0, c1)
1041 double precision, intent(in) :: c0, c1
1042 double precision, intent(inout) :: e_gas
1043 double precision :: xval, yval, der, deltax
1044 integer :: ii
1045
1046 yval = bigdouble
1047 xval = e_gas
1048 der = one
1049 deltax = one
1050 ii = 0
1051 !> Compare error with dx = dx/dy dy
1052 do while(dabs(deltax) .gt. fld_bisect_tol)
1053 yval = polynomial_bisection(xval, c0, c1)
1054 der = dpolynomial_bisection(xval, c0, c1)
1055 deltax = -yval/der
1056 xval = xval + deltax
1057 ii = ii + 1
1058 if(ii .gt. 1d3) then
1059 if(fld_debug)print*, 'skip to bisection algorithm'
1060 call bisection_method(e_gas, c0, c1)
1061 return
1062 endif
1063 enddo
1064 e_gas = xval
1065 end subroutine newton_method
1066
1067 !> Find the root of the 4th degree polynomial using the Halley method
1068 subroutine halley_method(e_gas, c0, c1)
1070 double precision, intent(in) :: c0, c1
1071 double precision, intent(inout) :: e_gas
1072 double precision :: xval, yval, der, dder, deltax
1073 integer :: ii
1074
1075 yval = bigdouble
1076 xval = e_gas
1077 der = one
1078 dder = one
1079 deltax = one
1080 ii = 0
1081 !> Compare error with dx = dx/dy dy
1082 do while (dabs(deltax) .gt. fld_bisect_tol)
1083 yval = polynomial_bisection(xval, c0, c1)
1084 der = dpolynomial_bisection(xval, c0, c1)
1085 dder = ddpolynomial_bisection(xval, c0, c1)
1086 deltax = -two*yval*der/(two*der**2 - yval*dder)
1087 xval = xval + deltax
1088 ii = ii + 1
1089 if(ii .gt. 1d3) then
1090 if(fld_debug)print*, 'skip to Newton algorithm'
1091 call newton_method(e_gas, c0, c1)
1092 return
1093 endif
1094 enddo
1095 e_gas = xval
1096 end subroutine halley_method
1097
1098 !> Evaluate polynomial at argument e_gas
1099 function polynomial_bisection(e_gas, c0, c1) result(val)
1101 double precision, intent(in) :: e_gas
1102 double precision, intent(in) :: c0, c1
1103 double precision :: val
1104
1105 val = e_gas**4.d0 + c1*e_gas - c0
1106 end function polynomial_bisection
1107
1108 !> Evaluate first derivative of polynomial at argument e_gas
1109 function dpolynomial_bisection(e_gas, c0, c1) result(der)
1111 double precision, intent(in) :: e_gas
1112 double precision, intent(in) :: c0, c1
1113 double precision :: der
1114
1115 der = 4.d0*e_gas**3.d0 + c1
1116 end function dpolynomial_bisection
1117
1118 !> Evaluate second derivative of polynomial at argument e_gas
1119 function ddpolynomial_bisection(e_gas, c0, c1) result(dder)
1121 double precision, intent(in) :: e_gas
1122 double precision, intent(in) :: c0, c1
1123 double precision :: dder
1124
1125 dder = 4.d0*3.d0*e_gas**2.d0
1126 end function ddpolynomial_bisection
1127end module mod_fld
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter one
Module for flux limited diffusion (FLD)-approximation in Radiation-(Magneto)hydrodynamics simulations...
Definition mod_fld.t:13
logical fld_no_mg
Definition mod_fld.t:33
subroutine, public fld_get_radpress(w, x, ixil, ixol, rad_pressure)
Returns Radiation Pressure as tensor NOTE: w is primitive on entry.
Definition mod_fld.t:476
double precision, public fld_bisect_tol
Tolerance for bisection method for Energy sourceterms This is a percentage of the minimum of gas- and...
Definition mod_fld.t:25
logical fld_force_mg_converged
switches for local debug purposes
Definition mod_fld.t:32
double precision, public fld_diff_tol
Tolerance for radiative Energy diffusion.
Definition mod_fld.t:27
double precision, public fld_gamma
A copy of (m)hd_gamma.
Definition mod_fld.t:46
subroutine fld_params_read(files)
public methods these are called in mod_hd_phys or mod_mhd_phys
Definition mod_fld.t:64
character(len=40) fld_fluxlimiter
flux limiter choice
Definition mod_fld.t:38
character(len=40) fld_opal_table
Definition mod_fld.t:36
subroutine update_diffcoeff(psa)
Definition mod_fld.t:858
double precision, public fld_boost_dt
Definition mod_fld.t:30
subroutine get_and_check_egas_erad_from_conserved(w, ixil, ixol, e_gas, e_rad)
Definition mod_fld.t:692
logical fld_slowsteps
handling ramp-up phase with explicit diffusion dt limit, slowly boosted
Definition mod_fld.t:29
subroutine, public fld_get_diffcoef_central(w, x, ixil, ixol)
Calculates cell-centered diffusion coefficient to be used in multigrid.
Definition mod_fld.t:932
subroutine evaluate_diffterm_onegrid(ixil, ixol, w, x)
inplace update of psa==>F_im(psa)
Definition mod_fld.t:895
double precision function ddpolynomial_bisection(e_gas, c0, c1)
Evaluate second derivative of polynomial at argument e_gas.
Definition mod_fld.t:1120
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...
Definition mod_fld.t:613
double precision, public fld_kappa0
Opacity value when using constant opacity.
Definition mod_fld.t:22
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 ...
Definition mod_fld.t:509
subroutine newton_method(e_gas, c0, c1)
Find the root of the 4th degree polynomial using the Newton method.
Definition mod_fld.t:1040
double precision function polynomial_bisection(e_gas, c0, c1)
Evaluate polynomial at argument e_gas.
Definition mod_fld.t:1100
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.
Definition mod_fld.t:644
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...
Definition mod_fld.t:437
character(len=40) fld_opacity_law
switches for opacity
Definition mod_fld.t:35
character(len=40) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
Definition mod_fld.t:44
subroutine fld_set_mg_bounds
Set the boundaries for the diffusion of E.
Definition mod_fld.t:150
subroutine fld_evaluate_implicit(qtc, psa)
inplace update of psa==>F_im(psa)
Definition mod_fld.t:875
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...
Definition mod_fld.t:197
logical fld_debug
Definition mod_fld.t:33
logical fld_radforce_split
source split for energy interact and radforce:
Definition mod_fld.t:18
double precision function dpolynomial_bisection(e_gas, c0, c1)
Evaluate first derivative of polynomial at argument e_gas.
Definition mod_fld.t:1110
subroutine, public fld_init(r_gamma)
Initialising FLD-module Read opacities Initialise Multigrid and adimensionalise kappa.
Definition mod_fld.t:86
logical fld_tiring_explicit
switch to handle photon tiring explicit or implicit
Definition mod_fld.t:20
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 ...
Definition mod_fld.t:531
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
Definition mod_fld.t:334
integer i_diff_mg
diffusion coefficient for multigrid method
Definition mod_fld.t:40
subroutine bisection_method(e_gas, c0, c1)
Find the root of the 4th degree polynomial using the bisection method.
Definition mod_fld.t:979
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...
Definition mod_fld.t:731
subroutine halley_method(e_gas, c0, c1)
Find the root of the 4th degree polynomial using the Halley method.
Definition mod_fld.t:1069
integer nth_for_diff_mg
diffusion coefficient stencil control
Definition mod_fld.t:42
Module with basic grid data structures.
Definition mod_forest.t:2
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
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 arad_norm
Normalised radiation constant.
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.
integer it_init
initial iteration count
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 slowsteps
If > 1, then in the first slowsteps-1 time steps dt is reduced by a factor .
integer, parameter bc_periodic
integer, parameter bc_special
boundary condition types
double precision c_norm
Normalised speed of light.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision unit_temperature
Physical scaling factor for temperature.
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
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)
integer refine_max_level
Maximal number of AMR levels.
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
Module for reading input and writing output.
subroutine get_global_minima(wmin, psa)
Compute global minima of iw variables over the leaves of the grid.
subroutine get_global_maxima(wmax, psa)
Compute global maxima of iw variables over the leaves of the grid.
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...
Definition mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:54
procedure(sub_get_tgas), pointer phys_get_tgas
Definition mod_physics.t:71
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
Definition mod_physics.t:81
procedure(sub_get_rfactor), pointer phys_get_rfactor
Definition mod_physics.t:73
integer phys_wider_stencil
To use wider stencils in flux calculations. A value of 1 will extend it by one cell in both direction...
Definition mod_physics.t:17
procedure(sub_set_mg_bounds), pointer phys_set_mg_bounds
Definition mod_physics.t:52
procedure(sub_implicit_update), pointer phys_implicit_update
Definition mod_physics.t:80
procedure(sub_get_csrad2), pointer phys_get_csrad2
Definition mod_physics.t:72
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...