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 !> Opacity per unit of unit_density
20 double precision, public :: fld_kappa0 = 0.34d0
21 !> Tolerance for bisection method for Energy sourceterms
22 !> This is a percentage of the minimum of gas- and radiation energy
23 double precision, public :: fld_bisect_tol = 1.d-4
24 !> Tolerance for radiative Energy diffusion
25 double precision, public :: fld_diff_tol = 1.d-4
26 !> switches for local debug purposes
28 !> switches for opacity
29 character(len=8) :: fld_opacity_law = 'const'
30 character(len=50) :: fld_opal_table = 'Y09800'
31 !> flux limiter choice
32 character(len=16) :: fld_fluxlimiter = 'Pomraning'
33 !> diffusion coefficient for multigrid method
34 integer :: i_diff_mg
35 !> diffusion coefficient stencil control
36 integer :: nth_for_diff_mg
37 !> Which method to find the root for the energy interaction polynomial
38 character(len=8) :: fld_interaction_method = 'Halley'
39 !> A copy of (m)hd_gamma
40 double precision, private, protected :: fld_gamma
41
42 !> public methods
43 !> these are called in mod_hd_phys or mod_mhd_phys
44 public :: fld_init
45 public :: fld_get_radpress
47 public :: add_fld_rad_force
48 public :: fld_radforce_get_dt
49 ! these are made public for mod_usr purposes and diagnostics
50 public :: fld_get_radflux
51 public :: fld_get_fluxlimiter
53 public :: fld_get_opacity_prim
54 contains
55
56 !> Reading in fld-list parameters from .par file
57 subroutine fld_params_read(files)
60 character(len=*), intent(in) :: files(:)
61 integer :: n
62
63 namelist /fld_list/ fld_kappa0, fld_radforce_split, &
66
67 do n = 1, size(files)
68 open(unitpar, file=trim(files(n)), status="old")
69 read(unitpar, fld_list, end=111)
70 111 close(unitpar)
71 end do
72 end subroutine fld_params_read
73
74 !> Initialising FLD-module
75 !> Read opacities
76 !> Initialise Multigrid and adimensionalise kappa
77 subroutine fld_init(He_abundance, r_gamma)
80 use mod_physics
83
84 double precision, intent(in) :: he_abundance, r_gamma
85 double precision :: sigma_thomson
86
88 fld_debug=.false.
89 fld_no_mg=.false.
91 if(.not.fld_no_mg)then
92 select case(nth_for_diff_mg)
93 case(1)
94 ! no need for stencil extension
95 case(2)
96 ! need for stencil extension
98 case default
99 call mpistop("nth_for_diff_mg must be 1 or 2")
100 end select
104 ! store the diffusion coefficient as extra variable (needed in mg vhelmholtz)
105 i_diff_mg = var_set_extravar("D", "D")
106 use_multigrid = .true.
107 ! use multigrid to solve a helmholtz equation with variable coefficient
108 mg%n_extra_vars = 1
109 mg%operator_type = mg_vhelmholtz
110 mg%smoother_type = mg_smoother_gs
111 endif
112 !> set gamma
113 fld_gamma = r_gamma
114 !> Read in opacity table if necesary
116 if(fld_opacity_law .eq. 'thomson') then
117 ! special fixed value for thomson opacity, assuming cgs units
118 sigma_thomson = 6.6524585d-25
119 fld_kappa0 = sigma_thomson/const_mp * (1.+2.*he_abundance)/(1.+4.*he_abundance)
120 endif
121 end subroutine fld_init
122
123 !> Set the boundaries for the diffusion of E
128
129 integer :: iB
130
131 ! Set boundary conditions for the multigrid solver
132 do ib = 1, 2*ndim
133 select case (typeboundary(iw_r_e, ib))
134 case (bc_symm)
135 ! d/dx u = 0
136 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
137 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
138 case (bc_asymm)
139 ! u = 0
140 mg%bc(ib, mg_iphi)%bc_type = mg_bc_dirichlet
141 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
142 case (bc_cont)
143 ! d/dx u = 0
144 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
145 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
146 case (bc_periodic)
147 ! Nothing to do here
148 case (bc_noinflow)
149 if (.not. associated(usr_special_mg_bc)) then
150 call mpistop("special BC for MG not defined")
151 else
152 call usr_special_mg_bc(ib)
153 endif
154 case (bc_special)
155 if (.not. associated(usr_special_mg_bc)) then
156 call mpistop("special BC for MG not defined")
157 else
158 call usr_special_mg_bc(ib)
159 endif
160 case default
161 call mpistop("divE_multigrid warning: unknown b.c. ")
162 end select
163 ! Neumann on the diffusion coefficient
164 ! d/dx u = 0
165 mg%bc(ib, mg_iveps)%bc_type = mg_bc_neumann
166 mg%bc(ib, mg_iveps)%bc_value = 0.0_dp
167 end do
168
169 end subroutine fld_set_mg_bounds
170
171
172
173 !> w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
174 !> This subroutine handles the radiation force and its work added explicitly
175 !> and the energy interaction term combined with photon tiring using an implicit update
176 subroutine add_fld_rad_force(qdt,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
177 use mod_constants
179 use mod_geometry
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
187
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
193
194 !> Calculate and add sourceterms
195 if(qsourcesplit .eqv. fld_radforce_split) then
196 active = .true.
197 nth_for_fld=2
198 ! store here lambda in a1 and fld_R in a2
199 call fld_get_eddington(wctprim,x,ixi^l,ixo^l,edd,a1,a2,nth_for_fld)
200 do idir = 1,ndim
201 call gradient(wctprim(ixi^s,iw_r_e),ixi^l,ixo^l,idir,tmp,nth_for_fld)
202 ! Radiation force = kappa*rho/c *Flux = lambda gradE
203 ! recycle grad E to store -lambda (grad E)_i
204 tmp(ixo^s) = -a1(ixo^s)*tmp(ixo^s)
205 !> Momentum equation source term
206 w(ixo^s,iw_mom(idir)) = w(ixo^s,iw_mom(idir))+ qdt*tmp(ixo^s)
207 !> Energy equation source term
208 w(ixo^s,iw_e) = w(ixo^s,iw_e) + qdt*wctprim(ixo^s,iw_mom(idir))*tmp(ixo^s)
209 enddo
210
211 !> Photon tiring : calculate tensor grad v (named div_v here)
212 ! NOTE: This is ok for uniform Cartesian only!!!!!
213 do idir = 1,ndim
214 do jdir = 1,ndim
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)
217 enddo
218 enddo
219 ! perform contraction fe : grad(v) with fe eddington tensor
220 {^ifoned
221 a3(ixo^s) = div_v(ixo^s,1,1)*edd(ixo^s,1,1)
222 }
223 {^iftwod
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)
228 }
229 {^ifthreed
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)
239 }
240 call fld_get_opacity_prim(wctprim,x,ixi^l,ixo^l,kappa)
241
243 sigma_b = const_rad_a*cc/4.d0*(unit_temperature**4.d0)/(unit_pressure)
244 !> e_gas is the INTERNAL ENERGY without KINETIC ENERGY
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)
248 endif
249 e_rad(ixo^s) = w(ixo^s,iw_r_e)
250
251 if(fix_small_values)then
252 {do ix^d = ixomin^d,ixomax^d\ }
253 e_gas(ix^d) = max(e_gas(ix^d),small_e)
254 e_rad(ix^d) = max(e_rad(ix^d),small_r_e)
255 {enddo\}
256 endif
257
258 !> Coefficients for the polynomial in Moens et al. 2022, eq 37. but with photon tiring (a3)
259 call phys_get_rfactor(wct,x,ixi^l,ixo^l,tmp)
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
263
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))
266
267 !> Loop over every cell for rootfinding method
268 {do ix^d = ixomin^d,ixomax^d\}
269 select case(fld_interaction_method)
270 case('Bisect')
271 call bisection_method(e_gas(ix^d),e_rad(ix^d),c0(ix^d),c1(ix^d))
272 case('Newton')
273 call newton_method(e_gas(ix^d),e_rad(ix^d),c0(ix^d),c1(ix^d))
274 case('Halley')
275 call halley_method(e_gas(ix^d),e_rad(ix^d),c0(ix^d),c1(ix^d))
276 case default
277 call mpistop('root-method not known')
278 end select
279 {enddo\}
280
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))
282
283 if(fix_small_values)then
284 {do ix^d = ixomin^d,ixomax^d\ }
285 e_gas(ix^d) = max(e_gas(ix^d),small_e)
286 e_rad(ix^d) = max(e_rad(ix^d),small_r_e)
287 {enddo\}
288 endif
289
290 !> Update gas-energy in w, internal + kinetic
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)
295 endif
296 !> Update rad-energy in w
297 w(ixo^s,iw_r_e) = e_rad(ixo^s)
298 end if
299 end subroutine add_fld_rad_force
300
301 !> get dt limit for radiation force and FLD explicit source additions
302 !> NOTE: w is primitive on entry
303 subroutine fld_radforce_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
306 use mod_geometry
308
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
312
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
318
319 if(fld_debug)print *,'DT limit on entry to radforce_get_dt=',dtnew
320 nth_for_fld=2
321 call fld_get_fluxlimiter_prim(w,x,ixi^l,ixo^l,lambda,fld_r,nth_for_fld)
322 if(slab_uniform) then
323 ^d&dxinv(^d)=one/dx^d;
324 do idim = 1, ndim
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)))
329 end do
330 else
331 do idim = 1, ndim
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))
336 end do
337 endif
338
339 if(fld_debug)print *,'DT limit after RADFORCE eff grav=',dtnew
340
341 ! here we interface back to fld_get_radpress
342 call phys_get_csrad2(w,x,ixi^l,ixo^l,tmp)
343 if(slab_uniform) then
344 ^d&dxinv(^d)=one/dx^d;
345 do idims=1,ndim
346 cmax(ixo^s)=dabs(w(ixo^s,iw_mom(idims)))+dsqrt(tmp(ixo^s))
347 if(idims==1) then
348 cmaxtot(ixo^s)=cmax(ixo^s)*dxinv(idims)
349 else
350 cmaxtot(ixo^s)=cmaxtot(ixo^s)+cmax(ixo^s)*dxinv(idims)
351 end if
352 end do
353 else
354 do idims=1,ndim
355 cmax(ixo^s)=dabs(w(ixo^s,iw_mom(idims)))+dsqrt(tmp(ixo^s))
356 if(idims==1) then
357 cmaxtot(ixo^s)=cmax(ixo^s)/block%ds(ixo^s,idims)
358 else
359 cmaxtot(ixo^s)=cmaxtot(ixo^s)+cmax(ixo^s)/block%ds(ixo^s,idims)
360 end if
361 end do
362 end if
363 ! courantmaxtots='max(summed c/dx)'
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
367
368 end subroutine fld_radforce_get_dt
369
370 !> Sets the opacity in the w-array
371 !> by calling mod_opal_opacity
372 !> NOTE: assumes primitives in w
373 !> NOTE: assuming opacity is local, not ok with cak line force
374 subroutine fld_get_opacity_prim(w, x, ixI^L, ixO^L, fld_kappa)
376 use mod_physics, only: phys_get_tgas
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)
383
384 integer :: ix^d
385 double precision :: rho0,temp0,n
386 double precision :: temp(ixi^s)
387
388 select case (fld_opacity_law)
389 case('const','thomson')
390 fld_kappa(ixo^s) = fld_kappa0/unit_opacity
391 case('opal')
392 call phys_get_tgas(w,x,ixi^l,ixo^l,temp)
393 {do ix^d=ixomin^d,ixomax^d\ }
394 rho0 = w(ix^d,iw_rho)*unit_density
395 temp0 = temp(ix^d)*unit_temperature
396 call set_opal_opacity(rho0,temp0,n)
397 fld_kappa(ix^d) = n/unit_opacity
398 {enddo\ }
399 case('special')
400 if (.not. associated(usr_special_opacity)) then
401 call mpistop("special opacity not defined")
402 endif
403 call usr_special_opacity(ixi^l, ixo^l, w, x, fld_kappa)
404 case default
405 call mpistop("Doesn't know opacity law")
406 end select
407 end subroutine fld_get_opacity_prim
408
409 !> Returns Radiation Pressure as tensor
410 !> NOTE: w is primitive on entry
411 subroutine fld_get_radpress(w, x, ixI^L, ixO^L, rad_pressure)
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)
417 integer :: i,j,nth
418 double precision :: eddington_tensor(ixi^s,1:ndim,1:ndim)
419 double precision :: lambda(ixi^s),fld_r(ixi^s)
420
421 nth=2
422 call fld_get_eddington(w, x, ixi^l, ixo^l, eddington_tensor, lambda, fld_r, nth)
423 if(fld_debug)then
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'
432 endif
433 do i=1,ndim
434 do j=1,ndim
435 rad_pressure(ixo^s,i,j) = eddington_tensor(ixo^s,i,j)*w(ixo^s,iw_r_e)
436 enddo
437 enddo
438 end subroutine fld_get_radpress
439
440 !> This subroutine calculates flux limiter lambda according to fld_fluxlimiter
441 !> It also calculates fld_R which is ratio of radiation scaleheight and mean free path
442 !> NOTE: nth and ixI and ixO not free to choose here: TODO
443 subroutine fld_get_fluxlimiter(w,x,ixI^L,ixO^L,fld_lambda,fld_R,nth)
445 use mod_geometry
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)
452
453 double precision :: wprim(ixi^s,1:nw)
454
455 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
456 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
457 call fld_get_fluxlimiter_prim(wprim,x,ixi^l,ixo^l,fld_lambda,fld_r,nth)
458
459 end subroutine fld_get_fluxlimiter
460
461 !> This subroutine calculates flux limiter lambda according to fld_fluxlimiter
462 !> It also calculates fld_R which is ratio of radiation scaleheight and mean free path
463 !> NOTE: this one operates on primitives
464 !> NOTE: nth and ixI and ixO not free to choose
465 subroutine fld_get_fluxlimiter_prim(w,x,ixI^L,ixO^L,fld_lambda,fld_R,nth)
467 use mod_geometry
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)
473
474 integer :: idir, ix^d
475 double precision :: kappa(ixi^s),normgrad2(ixi^s)
476 double precision :: grad_r_e(ixi^s)
477
478 select case(fld_fluxlimiter)
479 case('Diffusion')
480 ! optically thick limit
481 fld_lambda(ixo^s) = 1.d0/3.d0
482 fld_r(ixo^s) = zero
483 case('FreeStream')
484 ! optically thin limit
485 normgrad2(ixo^s) = zero
486 do idir=1,ndim
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
489 end do
490 call fld_get_opacity_prim(w,x,ixi^l,ixo^l,kappa)
491 ! Calculate R everywhere
492 ! |grad E|/(rho kappa E)
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)
495 ! treat uniform case as diffusion limit
496 fld_r(ixo^s)=zero
497 fld_lambda(ixo^s) = 1.0d0/3.0d0
498 elsewhere
499 fld_lambda(ixo^s) = one/fld_r(ixo^s)
500 endwhere
501 case('Pomraning')
502 ! Calculate R everywhere
503 ! |grad E|/(rho kappa E)
504 normgrad2(ixo^s) = zero
505 do idir = 1,ndim
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
508 end do
509 call fld_get_opacity_prim(w,x,ixi^l,ixo^l,kappa)
510 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
511 ! Calculate the flux limiter, lambda
512 ! Levermore and Pomraning: lambda = (2 + R)/(6 + 3R + R^2)
513 fld_lambda(ixo^s) = (2.d0+fld_r(ixo^s))/(6.d0+3*fld_r(ixo^s)+fld_r(ixo^s)**2)
514 case('Minerbo')
515 ! Calculate R everywhere
516 ! |grad E|/(rho kappa E)
517 normgrad2(ixo^s) = zero
518 do idir = 1,ndim
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
521 end do
522 call fld_get_opacity_prim(w, x, ixi^l, ixo^l, kappa)
523 fld_r(ixo^s) = dsqrt(normgrad2(ixo^s))/(kappa(ixo^s)*w(ixo^s,iw_rho)*w(ixo^s,iw_r_e))
524 ! Calculate the flux limiter, lambda
525 ! Minerbo:
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))
529 else
530 fld_lambda(ix^d) = 1.d0/(1.d0+fld_r(ix^d)+dsqrt(1.d0+2.d0*fld_r(ix^d)))
531 endif
532 {enddo\}
533 case('special')
534 if (.not. associated(usr_special_fluxlimiter)) then
535 call mpistop("special fluxlimiter not defined")
536 endif
537 call usr_special_fluxlimiter(ixi^l,ixo^l,w,x,fld_lambda,fld_r)
538 case default
539 call mpistop('Fluxlimiter unknown')
540 end select
541
542 end subroutine fld_get_fluxlimiter_prim
543
544 !> Calculate Radiation Flux (use of cgs units)
545 !> NOTE: only for diagnostics purposes (w conservative on entry)
546 !> This returns cell centered values for radiation flux
547 subroutine fld_get_radflux(w, x, ixI^L, ixO^L, rad_flux)
549 use mod_geometry
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)
555
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)
561
562 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
563 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
564
565 cc = const_c/unit_velocity
566 call fld_get_opacity_prim(wprim, x, ixi^l, ixo^l, kappa)
567 nth_for_fld=2
568 call fld_get_fluxlimiter_prim(wprim, x, ixi^l, ixo^l, lambda, fld_r, nth_for_fld)
569 !> Calculate the Flux using the fld closure relation
570 !> F = -c*lambda/(kappa*rho) *grad E
571 do idir = 1,ndim
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)
574 end do
575 end subroutine fld_get_radflux
576
577 !> Calculate Eddington-tensor (where w is primitive)
578 !> also feeds back the flux limiter lambda and R
579 subroutine fld_get_eddington(w, x, ixI^L, ixO^L, eddington_tensor, lambda, fld_R, nth)
581 use mod_geometry
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)
587
588 integer :: idir,jdir
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)
592
593 normgrad2(ixo^s) = zero
594 do idir = 1,ndim
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
598 end do
599 do idir = 1,ndim
600 do jdir = 1,ndim
601 if(idir==jdir)then
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)
603 else
604 nn_regularized(ixo^s,idir,jdir)=(grad_r_e(ixo^s,idir)*grad_r_e(ixo^s,jdir))/(normgrad2(ixo^s)+smalldouble**2)
605 endif
606 enddo
607 enddo
608 ! get lambda and R
609 call fld_get_fluxlimiter_prim(w,x,ixi^l,ixo^l,lambda,fld_r,nth)
610 ! store f_e= lambda + lambda^2 R^2
611 tmp(ixo^s) = lambda(ixo^s)+(lambda(ixo^s)*fld_r(ixo^s))**2
612 do idir = 1,ndim
613 ! first compute the isotropic (diagonal) part
614 eddington_tensor(ixo^s,idir,idir) = half*(one-tmp(ixo^s))
615 enddo
616 do idir = 1,ndim
617 do jdir = 1,ndim
618 ! initialize off-diagonal part here
619 if(idir .ne. jdir) eddington_tensor(ixo^s,idir,jdir) = zero
620 ! add part depending on unit vectors along gradient E
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)
623 enddo
624 enddo
625 end subroutine fld_get_eddington
626
627 !> Calling all subroutines to perform the multigrid method
628 !> Communicates rad_e and diff_coeff to multigrid library
629 !> Advance psa=psb+dtfactor*qdt*F_im(psa)
630 subroutine fld_implicit_update(dtfactor,qdt,qtC,psa,psb)
632 use mod_forest
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
639
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
644
645 ! Avoid setting a very restrictive limit to residual when time step small
646 if(qdt < dtmin) then
647 if(mype==0)then
648 print *,'skipping implicit solve: dt too small!'
649 print *,'Currently at time=',global_time,' time step=',qdt,' dtmin=',dtmin
650 endif
651 return
652 endif
653 max_residual = fld_diff_tol
654
655 call mg_set_methods(mg)
656 if(.not. mg%is_allocated) call mpistop("multigrid tree not allocated yet")
657
658 ! Here we handle the global helmholtz problem with variable coefficient
659 ! The equation we solve is div(D nabla Erad^(n+1)) -(1/dt)Erad^(n+1)=-(1/dt)Erad^n
660 ! Helmholtz equation is div(eps nabla phi) -lambda phi = f
661 ! hence eps is our variable coefficient D=fld_lambda*c/(kappa*rho)
662 ! hence phi is Erad and lambda=1/dt
663 lambda = 1.d0/(dtfactor *qdt)
664 fac = 1.d0
665 call vhelmholtz_set_lambda(lambda)
666 call update_diffcoeff(psa)
667 ! copy in the (variable) diffusion coefficient in mg_iveps
668 call mg_copy_to_tree(i_diff_mg, mg_iveps, factor=fac, state_from=psa)
669 ! copy in the Erad variable in mg_iphi (the one we solve for)
670 call mg_copy_to_tree(iw_r_e, mg_iphi, factor=fac, state_from=psa)
671 ! copy in RHS f factor as Erad with factor -(1/dt)
672 call mg_copy_to_tree(iw_r_e, mg_irhs, factor=-lambda, state_from=psb)
673 if(time_advance)then
674 call mg_restrict(mg, mg_iveps)
675 call mg_fill_ghost_cells(mg, mg_iveps)
676 endif
677 call mg_fas_fmg(mg, .true., max_res=res)
678 do n = 1, max_its
679 if(res < max_residual) exit
680 call mg_fas_vcycle(mg, max_res=res)
681 end do
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
688 endif
689 call mpistop("no convergence in MG")
690 end if
691 if(res .eq. 0.d0) then
692 if (mg%my_rank == 0) write(*,*) it, 'ZERO residual from MG ', res
693 endif
694 ! copy back the Erad variable in iw_r_e
695 call mg_copy_from_tree_gc(mg_iphi, iw_r_e, state_to=psa)
696
697 if(fix_small_values)then
698 ixo^l=ixm^ll;
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)
702 {enddo\}
703 end do
704 endif
705
706 end subroutine fld_implicit_update
707
708 subroutine update_diffcoeff(psa)
710 type(state), target :: psa(max_blocks)
711
712 integer :: iigrid, igrid, ixO^L
713
714 ! we will need diffusion coefficients in i-1 i i+1
715 ixo^l=ixm^ll^ladd1;
716 !$OMP PARALLEL DO PRIVATE(igrid)
717 do iigrid=1,igridstail; igrid=igrids(iigrid);
718 call fld_get_diffcoef_central(psa(igrid)%w, psa(igrid)%x, ixg^ll, ixo^l)
719 end do
720 !$OMP END PARALLEL DO
721
722 end subroutine update_diffcoeff
723
724 !> inplace update of psa==>F_im(psa)
725 subroutine fld_evaluate_implicit(qtC,psa)
727 type(state), target :: psa(max_blocks)
728 double precision, intent(in) :: qtC
729 integer :: iigrid, igrid
730 integer :: ixO^L
731
732 call update_diffcoeff(psa)
733
734 ixo^l=ixm^ll;
735 !$OMP PARALLEL DO PRIVATE(igrid)
736 do iigrid=1,igridstail; igrid=igrids(iigrid);
737 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
738 call evaluate_diffterm_onegrid(ixg^ll,ixo^l,psa(igrid)%w)
739 end do
740 !$OMP END PARALLEL DO
741
742 end subroutine fld_evaluate_implicit
743
744 !> inplace update of psa==>F_im(psa)
745 subroutine evaluate_diffterm_onegrid(ixI^L,ixO^L,w)
747 integer, intent(in) :: ixI^L, ixO^L
748 double precision, intent(inout) :: w(ixI^S, 1:nw)
749
750 double precision :: divF(ixI^S)
751 integer :: idir, jxO^L, hxO^L
752
753 if(.not.slab)call mpistop("laplacian coded up for uniform cartesian grid")
754
755 ! Here we use diffusion coefficients in positions i-1 i i+1 and exploit harmonic means i.e. 2ab/(a+b)
756 ! since this is how the multigrid library handles the div(eps nabla phi) term in Helmholtz equation
757 ! div(eps nabla phi) - lambda phi = f
758 divf(ixo^s) = 0.d0
759 do idir = 1,ndim
760 hxo^l=ixo^l-kr(idir,^d);
761 jxo^l=ixo^l+kr(idir,^d);
762 divf(ixo^s) = divf(ixo^s) + &
763 (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)) &
764 -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)) &
765 +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))) &
766 +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
767 ! below uses artihmetic mean, different from mg method
768 !divF(ixO^S) = divF(ixO^S) + &
769 ! (w(jxO^S,iw_r_e)*half*(w(ixO^S,i_diff_mg) + w(jxO^S,i_diff_mg)) &
770 ! -w(ixO^S,iw_r_e)*(half*(w(ixO^S,i_diff_mg) + w(jxO^S,i_diff_mg)) &
771 ! +half*(w(ixO^S,i_diff_mg) + w(hxO^S,i_diff_mg))) &
772 ! +w(hxO^S,iw_r_e)*half*(w(ixO^S,i_diff_mg) + w(hxO^S,i_diff_mg)))/dxlevel(idir)**2
773 enddo
774 w(ixo^s,iw_r_e) = divf(ixo^s)
775
776 end subroutine evaluate_diffterm_onegrid
777
778 !> Calculates cell-centered diffusion coefficient to be used in multigrid
779 subroutine fld_get_diffcoef_central(w, x, ixI^L, ixO^L)
781 use mod_geometry
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)
787
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)
791
792 cc = const_c/unit_velocity
793 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
794 ! ensure entries in entire ixI range
795 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
796 call fld_get_opacity_prim(wprim, x, ixi^l, ixo^l, kappa)
797 ! note that we use central difference here (last argument is 1 or 2)
798 call fld_get_fluxlimiter_prim(wprim, x, ixi^l, ixo^l, lambda, fld_r, nth_for_diff_mg)
799 w(ixo^s,i_diff_mg) = cc*lambda(ixo^s)/(kappa(ixo^s)*wprim(ixo^s,iw_rho))
800 if(associated(usr_special_diffcoef)) call usr_special_diffcoef(w, wprim, x, ixi^l, ixo^l)
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")
804 endif
805 if(maxval(w(ixo^s,i_diff_mg))>bigdouble) call mpistop("too large diffusion coefficient")
806
807 if(fld_debug)then
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
822 endif
823
824 end subroutine fld_get_diffcoef_central
825
826
827 !> Find the root of the 4th degree polynomial using the bisection method
828 subroutine bisection_method(e_gas, E_rad, c0, c1)
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
835
836 n = 0
837 max_its = 100
838 bisect_a = zero
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
842 n = n +1
843 if(n .gt. max_its) then
844 call mpistop('No convergece in bisection scheme')
845 endif
846 if(polynomial_bisection(bisect_a, c0, c1)*&
847 polynomial_bisection(bisect_b, c0, c1) .lt. zero) then
848 if(polynomial_bisection(bisect_a, c0, c1)*&
849 polynomial_bisection(bisect_c, c0, c1) .lt. zero) then
850 bisect_b = bisect_c
851 elseif(polynomial_bisection(bisect_b, c0, c1)*&
852 polynomial_bisection(bisect_c, c0, c1) .lt. zero) then
853 bisect_a = bisect_c
854 elseif(polynomial_bisection(bisect_a, c0, c1) .eq. zero) then
855 bisect_b = bisect_a
856 bisect_c = bisect_a
857 goto 2435
858 elseif(polynomial_bisection(bisect_b, c0, c1) .eq. zero) then
859 bisect_a = bisect_b
860 bisect_c = bisect_b
861 goto 2435
862 elseif(polynomial_bisection(bisect_c, c0, c1) .eq. zero) then
863 bisect_a = bisect_c
864 bisect_b = bisect_c
865 goto 2435
866 else
867 call mpistop("Problem with fld bisection method")
868 endif
869 elseif(polynomial_bisection(bisect_a, c0, c1) &
870 - polynomial_bisection(bisect_b, c0, c1) .lt. fld_bisect_tol*polynomial_bisection(bisect_a, c0, c1)) then
871 goto 2435
872 else
873 bisect_a = e_gas
874 bisect_b = e_gas
875 if(fld_debug)print*, "IGNORING GAS-RAD ENERGY EXCHANGE ", c0, c1
876 if(fld_debug)print*, polynomial_bisection(bisect_a, c0, c1), polynomial_bisection(bisect_b, c0, c1)
877 call mpistop('issues in bisection scheme')
878 if(polynomial_bisection(bisect_a, c0, c1) .le. smalldouble) then
879 bisect_b = bisect_a
880 elseif(polynomial_bisection(bisect_a, c0, c1) .le. smalldouble) then
881 bisect_a = bisect_b
882 endif
883 goto 2435
884 endif
885 enddo
886 2435 e_gas = (bisect_a + bisect_b)/two
887 end subroutine bisection_method
888
889 !> Find the root of the 4th degree polynomial using the Newton method
890 subroutine newton_method(e_gas, E_rad, c0, c1)
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
896 integer :: ii
897
898 yval = bigdouble
899 xval = e_gas
900 der = one
901 deltax = one
902 ii = 0
903 !> Compare error with dx = dx/dy dy
904 do while(dabs(deltax) .gt. fld_bisect_tol)
905 yval = polynomial_bisection(xval, c0, c1)
906 der = dpolynomial_bisection(xval, c0, c1)
907 deltax = -yval/der
908 xval = xval + deltax
909 ii = ii + 1
910 if(ii .gt. 1d3) then
911 if(fld_debug)print*, 'skip to bisection algorithm'
912 call bisection_method(e_gas, e_rad, c0, c1)
913 return
914 endif
915 enddo
916 e_gas = xval
917 end subroutine newton_method
918
919 !> Find the root of the 4th degree polynomial using the Halley method
920 subroutine halley_method(e_gas, E_rad, c0, c1)
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
926 integer :: ii
927
928 yval = bigdouble
929 xval = e_gas
930 der = one
931 dder = one
932 deltax = one
933 ii = 0
934 !> Compare error with dx = dx/dy dy
935 do while (dabs(deltax) .gt. fld_bisect_tol)
936 yval = polynomial_bisection(xval, c0, c1)
937 der = dpolynomial_bisection(xval, c0, c1)
938 dder = ddpolynomial_bisection(xval, c0, c1)
939 deltax = -two*yval*der/(two*der**2 - yval*dder)
940 xval = xval + deltax
941 ii = ii + 1
942 if(ii .gt. 1d3) then
943 if(fld_debug)print*, 'skip to Newton algorithm'
944 !call mpistop("iterates to 1000, entering Newton method")
945 call newton_method(e_gas, e_rad, c0, c1)
946 return
947 endif
948 enddo
949 e_gas = xval
950 end subroutine halley_method
951
952 !> Evaluate polynomial at argument e_gas
953 function polynomial_bisection(e_gas, c0, c1) result(val)
955 double precision, intent(in) :: e_gas
956 double precision, intent(in) :: c0, c1
957 double precision :: val
958
959 val = e_gas**4.d0 + c1*e_gas - c0
960 end function polynomial_bisection
961
962 !> Evaluate first derivative of polynomial at argument e_gas
963 function dpolynomial_bisection(e_gas, c0, c1) result(der)
965 double precision, intent(in) :: e_gas
966 double precision, intent(in) :: c0, c1
967 double precision :: der
968
969 der = 4.d0*e_gas**3.d0 + c1
970 end function dpolynomial_bisection
971
972 !> Evaluate second derivative of polynomial at argument e_gas
973 function ddpolynomial_bisection(e_gas, c0, c1) result(dder)
975 double precision, intent(in) :: e_gas
976 double precision, intent(in) :: c0, c1
977 double precision :: dder
978
979 dder = 4.d0*3.d0*e_gas**2.d0
980 end function ddpolynomial_bisection
981end module mod_fld
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...
Definition mod_fld.t:13
logical fld_no_mg
Definition mod_fld.t:27
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:412
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:23
double precision, public fld_diff_tol
Tolerance for radiative Energy diffusion.
Definition mod_fld.t:25
subroutine fld_params_read(files)
public methods these are called in mod_hd_phys or mod_mhd_phys
Definition mod_fld.t:58
subroutine update_diffcoeff(psa)
Definition mod_fld.t:709
character(len=8) fld_opacity_law
switches for opacity
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:780
subroutine evaluate_diffterm_onegrid(ixil, ixol, w)
inplace update of psa==>F_im(psa)
Definition mod_fld.t:746
double precision function ddpolynomial_bisection(e_gas, c0, c1)
Evaluate second derivative of polynomial at argument e_gas.
Definition mod_fld.t:974
character(len=16) fld_fluxlimiter
flux limiter choice
Definition mod_fld.t:32
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...
Definition mod_fld.t:548
subroutine bisection_method(e_gas, e_rad, c0, c1)
Find the root of the 4th degree polynomial using the bisection method.
Definition mod_fld.t:829
character(len=50) fld_opal_table
Definition mod_fld.t:30
double precision, public fld_kappa0
Opacity per unit of unit_density.
Definition mod_fld.t:20
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:444
double precision function polynomial_bisection(e_gas, c0, c1)
Evaluate polynomial at argument e_gas.
Definition mod_fld.t:954
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:580
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:375
subroutine, public fld_init(he_abundance, r_gamma)
Initialising FLD-module Read opacities Initialise Multigrid and adimensionalise kappa.
Definition mod_fld.t:78
subroutine fld_set_mg_bounds
Set the boundaries for the diffusion of E.
Definition mod_fld.t:125
subroutine halley_method(e_gas, e_rad, c0, c1)
Find the root of the 4th degree polynomial using the Halley method.
Definition mod_fld.t:921
subroutine fld_evaluate_implicit(qtc, psa)
inplace update of psa==>F_im(psa)
Definition mod_fld.t:726
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:177
character(len=8) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
Definition mod_fld.t:38
subroutine newton_method(e_gas, e_rad, c0, c1)
Find the root of the 4th degree polynomial using the Newton method.
Definition mod_fld.t:891
logical fld_debug
switches for local debug purposes
Definition mod_fld.t:27
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:964
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:466
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:304
integer i_diff_mg
diffusion coefficient for multigrid method
Definition mod_fld.t:34
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:631
integer nth_for_diff_mg
diffusion coefficient stencil control
Definition mod_fld.t:36
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.
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
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, dimension(:,:), allocatable rnode
Corner coordinates.
double precision unit_temperature
Physical scaling factor for temperature.
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)
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...
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...