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