35 real(8),
allocatable,
private :: ay(:), wy(:), aphi(:), wphi(:)
38 real(8),
private :: cak_gamma
41 real(8),
private :: lum, dlum, drstar, dke, dclight
44 real(8),
private :: tfloor
50 integer,
parameter,
private :: radstream=0, fdisc=1, fdisc_cutoff=2
79 character(len=*),
intent(in) :: files(:)
89 open(
unitpar, file=trim(files(n)), status=
"old")
90 read(
unitpar, cak_list,
end=111)
101 real(8),
intent(in) :: phys_gamma
103 cak_gamma = phys_gamma
116 gcak1_ = var_set_extravar(
"gcak1",
"gcak1")
117 fdf_ = var_set_extravar(
"fdfac",
"fdfac")
121 gcak1_ = var_set_extravar(
"gcak1",
"gcak1")
122 gcak2_ = var_set_extravar(
"gcak2",
"gcak2")
123 gcak3_ = var_set_extravar(
"gcak3",
"gcak3")
129 call mpistop(
'CAK error: choose alpha in [0,1[')
133 call mpistop(
'CAK error: chosen Qbar or Q0 is < 0')
137 call mpistop(
'CAK error: choose either 1-D or vector force')
147 real(8),
intent(in) :: rstar, twind
148 double precision :: const_kappae_local
149 const_kappae_local=0.34d0
165 integer,
intent(in) :: ixI^L, ixO^L
166 real(8),
intent(in) :: qdt, x(ixI^S,1:ndim), wCT(ixI^S,1:nw)
167 real(8),
intent(inout) :: w(ixI^S,1:nw)
168 logical,
intent(in) :: energy, qsourcesplit
169 logical,
intent(inout) :: active
172 real(8) :: gl(ixO^S,1:3), ge(ixO^S), ptherm(ixI^S), pmin(ixI^S)
184 gl(ixo^s,1:3) = 0.0d0
191 call mpistop(
"No valid force option")
196 if (idir == 1) gl(ixo^s,idir) = gl(ixo^s,idir) + ge(ixo^s)
198 w(ixo^s,iw_mom(idir)) = w(ixo^s,iw_mom(idir)) &
199 + qdt * gl(ixo^s,idir) * wct(ixo^s,iw_rho)
202 w(ixo^s,iw_e) = w(ixo^s,iw_e) + qdt * gl(ixo^s,idir) * wct(ixo^s,iw_mom(idir))
208 call phys_get_pthermal(w,x,ixi^l,ixo^l,ptherm)
209 pmin(ixo^s) = w(ixo^s,iw_rho) * tfloor
211 where (ptherm(ixo^s) < pmin(ixo^s))
212 w(ixo^s,iw_e) = w(ixo^s,iw_e) + (pmin(ixo^s) - ptherm(ixo^s))/(cak_gamma - 1.0d0)
224 integer,
intent(in) :: ixI^L, ixO^L
225 real(8),
intent(in) :: wCT(ixI^S,1:nw), x(ixI^S,1:ndim)
226 real(8),
intent(inout) :: w(ixI^S,1:nw)
227 real(8),
intent(inout) :: gcak(ixO^S,1:3)
230 real(8) :: vr(ixI^S), dvrdr(ixO^S)
231 real(8) :: beta_fd(ixO^S), fdfac(ixO^S), taus(ixO^S), ge(ixO^S)
233 vr(ixi^s) = wct(ixi^s,iw_mom(1)) / wct(ixi^s,iw_rho)
236 if (physics_type ==
'hd')
then
238 dvrdr(ixo^s) = abs(dvrdr(ixo^s))
241 dvrdr(ixo^s) = max(dvrdr(ixo^s), smalldouble)
249 case(radstream, fdisc)
250 taus(ixo^s) =
gayley_qbar * dke * dclight * wct(ixo^s,iw_rho)/dvrdr(ixo^s)
255 taus(ixo^s) =
gayley_q0 * dke * dclight * wct(ixo^s,iw_rho)/dvrdr(ixo^s)
257 * ( (1.0d0 + taus(ixo^s))**(1.0d0 -
cak_alpha) - 1.0d0 ) &
260 call mpistop(
"Error in force computation.")
264 beta_fd(ixo^s) = ( 1.0d0 - vr(ixo^s)/(x(ixo^s,1) * dvrdr(ixo^s)) ) &
265 * (drstar/x(ixo^s,1))**2.0d0
270 case(fdisc, fdisc_cutoff)
271 where (beta_fd(ixo^s) >= 1.0d0)
273 elsewhere (beta_fd(ixo^s) < -1.0d10)
275 elsewhere (abs(beta_fd(ixo^s)) > 1.0
d-3)
276 fdfac(ixo^s) = (1.0d0 - (1.0d0 - beta_fd(ixo^s))**(1.0d0 +
cak_alpha)) &
279 fdfac(ixo^s) = 1.0d0 - 0.5d0*
cak_alpha*beta_fd(ixo^s) &
280 * (1.0d0 + 1.0d0/3.0d0 * (1.0d0 -
cak_alpha)*beta_fd(ixo^s))
285 gcak(ixo^s,1) = gcak(ixo^s,1) * fdfac(ixo^s)
286 gcak(ixo^s,2) = 0.0d0
287 gcak(ixo^s,3) = 0.0d0
290 w(ixo^s,
gcak1_) = gcak(ixo^s,1)
291 w(ixo^s,
fdf_) = fdfac(ixo^s)
301 integer,
intent(in) :: ixI^L, ixO^L
302 real(8),
intent(in) :: wCT(ixI^S,1:nw), x(ixI^S,1:ndim)
303 real(8),
intent(inout) :: w(ixI^S,1:nw)
304 real(8),
intent(inout) :: gcak(ixO^S,1:3)
307 real(8) :: a1, a2, a3, wyray, y, wpray, phiray, wtot, mustar, dvndn
308 real(8) :: costp, costp2, sintp, cospp, sinpp, cott0
309 real(8) :: vr(ixI^S), vt(ixI^S), vp(ixI^S)
310 real(8) :: vrr(ixI^S), vtr(ixI^S), vpr(ixI^S)
311 real(8) :: dvrdr(ixO^S), dvtdr(ixO^S), dvpdr(ixO^S)
312 real(8) :: dvrdt(ixO^S), dvtdt(ixO^S), dvpdt(ixO^S)
313 real(8) :: dvrdp(ixO^S), dvtdp(ixO^S), dvpdp(ixO^S)
314 integer :: ix^D, itray, ipray
317 vt(ixo^s) = 0.0d0; vtr(ixo^s) = 0.0d0
318 vp(ixo^s) = 0.0d0; vpr(ixo^s) = 0.0d0
320 dvrdr(ixo^s) = 0.0d0; dvtdr(ixo^s) = 0.0d0; dvpdr(ixo^s) = 0.0d0
321 dvrdt(ixo^s) = 0.0d0; dvtdt(ixo^s) = 0.0d0; dvpdt(ixo^s) = 0.0d0
322 dvrdp(ixo^s) = 0.0d0; dvtdp(ixo^s) = 0.0d0; dvpdp(ixo^s) = 0.0d0
325 vr(ixi^s) = wct(ixi^s,iw_mom(1)) / wct(ixi^s,iw_rho)
326 vrr(ixi^s) = vr(ixi^s) / x(ixi^s,1)
329 vt(ixi^s) = wct(ixi^s,iw_mom(2)) / wct(ixi^s,iw_rho)
330 vtr(ixi^s) = vt(ixi^s) / x(ixi^s,1)
333 vp(ixi^s) = wct(ixi^s,iw_mom(3)) / wct(ixi^s,iw_rho)
334 vpr(ixi^s) = vp(ixi^s) / x(ixi^s,1)
358 {
do ix^db=ixomin^db,ixomax^db\}
379 mustar = sqrt(max(1.0d0 - (drstar/x(ix^d,1))**2.0d0, 0.0d0))
380 costp = 1.0d0 - y*(1.0d0 - mustar)
382 sintp = sqrt(max(1.0d0 - costp2, 0.0d0))
385 {^nooned cott0 = cos(x(ix^d,2))/sin(x(ix^d,2))}
388 wtot = wyray * wpray * (1.0d0 - mustar)
396 dvndn = a1*a1 * dvrdr(ix^d) + a2*a2 * (dvtdt(ix^d) + vrr(ix^d)) &
397 + a3*a3 * (dvpdp(ix^d) + cott0 * vtr(ix^d) + vrr(ix^d)) &
398 + a1*a2 * (dvtdr(ix^d) + dvrdt(ix^d) - vtr(ix^d)) &
399 + a1*a3 * (dvpdr(ix^d) + dvrdp(ix^d) - vpr(ix^d)) &
400 + a2*a3 * (dvpdt(ix^d) + dvtdp(ix^d) - cott0 * vpr(ix^d))
407 gcak(ix^d,1) = gcak(ix^d,1) + (dvndn/wct(ix^d,iw_rho))**
cak_alpha * a1 * wtot
408 gcak(ix^d,2) = gcak(ix^d,2) + (dvndn/wct(ix^d,iw_rho))**
cak_alpha * a2 * wtot
409 gcak(ix^d,3) = gcak(ix^d,3) + (dvndn/wct(ix^d,iw_rho))**
cak_alpha * a3 * wtot
417 * dlum/(4.0d0*dpi*drstar**2.0d0 * dclight**(1.0d0+
cak_alpha)) &
421 gcak(ixo^s,2) = 0.0d0
422 gcak(ixo^s,3) = 0.0d0
426 w(ixo^s,
gcak1_) = gcak(ixo^s,1)
427 w(ixo^s,
gcak2_) = gcak(ixo^s,2)
428 w(ixo^s,
gcak3_) = gcak(ixo^s,3)
436 integer,
intent(in) :: ixI^L, ixO^L
437 real(8),
intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
438 real(8),
intent(out):: ge(ixO^S)
440 ge(ixo^s) = dke * dlum/(4.0d0*dpi * dclight * x(ixo^s,1)**2.0d0)
448 integer,
intent(in) :: ixI^L, ixO^L
449 real(8),
intent(in) :: dx^D, x(ixI^S,1:ndim)
450 real(8),
intent(in) :: wprim(ixI^S,1:nw)
451 real(8),
intent(inout) :: dtnew
454 real(8) :: ge(ixO^S), max_gr, dt_cak
461 max_gr = max( maxval(abs(ge(ixo^s) + wprim(ixo^s,
gcak1_))), epsilon(1.0d0) )
462 dt_cak = minval( sqrt(
block%dx(ixo^s,1)/max_gr) )
467 max_gr = max( maxval(abs(wprim(ixo^s,
gcak2_))), epsilon(1.0d0) )
468 dt_cak = minval( sqrt(
block%dx(ixo^s,1) *
block%dx(ixo^s,2)/max_gr) )
472 max_gr = max( maxval(abs(wprim(ixo^s,
gcak3_))), epsilon(1.0d0) )
473 dt_cak = minval( sqrt(
block%dx(ixo^s,1) * sin(
block%dx(ixo^s,3))/max_gr) )
485 integer,
intent(in) :: ixI^L, ixO^L, idir
486 real(8),
intent(in) :: vfield(ixI^S), x(ixI^S,1:ndim)
487 real(8),
intent(out) :: grad_vn(ixO^S)
490 real(8) :: forw(ixO^S), backw(ixO^S), cent(ixO^S)
491 integer :: jrx^L, hrx^L{^NOONED,jtx^L, htx^L}{^IFTHREED,jpx^L, hpx^L}
494 jrx^l=ixo^l+
kr(1,^
d);
495 hrx^l=ixo^l-
kr(1,^
d);
499 jtx^l=ixo^l+
kr(2,^
d);
500 htx^l=ixo^l-
kr(2,^
d);
505 jpx^l=ixo^l+
kr(3,^
d);
506 hpx^l=ixo^l-
kr(3,^
d);
512 forw(ixo^s) = (x(ixo^s,1) - x(hrx^s,1)) * vfield(jrx^s) &
513 / ((x(jrx^s,1) - x(ixo^s,1)) * (x(jrx^s,1) - x(hrx^s,1)))
515 backw(ixo^s) = -(x(jrx^s,1) - x(ixo^s,1)) * vfield(hrx^s) &
516 / ((x(ixo^s,1) - x(hrx^s,1)) * (x(jrx^s,1) - x(hrx^s,1)))
518 cent(ixo^s) = (x(jrx^s,1) + x(hrx^s,1) - 2.0d0*x(ixo^s,1)) * vfield(ixo^s) &
519 / ((x(ixo^s,1) - x(hrx^s,1)) * (x(jrx^s,1) - x(ixo^s,1)))
522 forw(ixo^s) = (x(ixo^s,2) - x(htx^s,2)) * vfield(jtx^s) &
523 / (x(ixo^s,1) * (x(jtx^s,2) - x(ixo^s,2)) * (x(jtx^s,2) - x(htx^s,2)))
525 backw(ixo^s) = -(x(jtx^s,2) - x(ixo^s,2)) * vfield(htx^s) &
526 / ( x(ixo^s,1) * (x(ixo^s,2) - x(htx^s,2)) * (x(jtx^s,2) - x(htx^s,2)))
528 cent(ixo^s) = (x(jtx^s,2) + x(htx^s,2) - 2.0d0*x(ixo^s,2)) * vfield(ixo^s) &
529 / ( x(ixo^s,1) * (x(ixo^s,2) - x(htx^s,2)) * (x(jtx^s,2) - x(ixo^s,2)))
533 forw(ixo^s) = (x(ixo^s,3) - x(hpx^s,3)) * vfield(jpx^s) &
534 / ( x(ixo^s,1)*sin(x(ixo^s,2)) * (x(jpx^s,3) - x(ixo^s,3)) * (x(jpx^s,3) - x(hpx^s,3)))
536 backw(ixo^s) = -(x(jpx^s,3) - x(ixo^s,3)) * vfield(hpx^s) &
537 / ( x(ixo^s,1)*sin(x(ixo^s,2)) * (x(ixo^s,3) - x(hpx^s,3)) * (x(jpx^s,3) - x(hpx^s,3)))
539 cent(ixo^s) = (x(jpx^s,3) + x(hpx^s,3) - 2.0d0*x(ixo^s,3)) * vfield(ixo^s) &
540 / ( x(ixo^s,1)*sin(x(ixo^s,2)) * (x(ixo^s,3) - x(hpx^s,3)) * (x(jpx^s,3) - x(ixo^s,3)))
545 grad_vn(ixo^s) = backw(ixo^s) + cent(ixo^s) + forw(ixo^s)
554 integer,
intent(in) :: ntheta_point, nphi_point
557 real(8) :: ymin, ymax, phipmin, phipmax, adum
569 allocate(ay(ntheta_point))
570 allocate(wy(ntheta_point))
571 allocate(aphi(nphi_point))
572 allocate(wphi(nphi_point))
597 print*,
'==========================='
598 print*,
' Radiation ray setup '
599 print*,
'==========================='
600 print*,
'Theta ray points + weights '
601 do ii = 1,ntheta_point
602 print*,ii,ay(ii),wy(ii)
605 print*,
'Phi ray points + weights '
607 print*,ii,aphi(ii),wphi(ii)
618 call mpi_bcast(ntheta_point,1,mpi_integer,0,
icomm,
ierrmpi)
619 call mpi_bcast(nphi_point,1,mpi_integer,0,
icomm,
ierrmpi)
622 allocate(ay(ntheta_point))
623 allocate(wy(ntheta_point))
624 allocate(aphi(nphi_point))
625 allocate(wphi(nphi_point))
628 call mpi_bcast(ay,ntheta_point,mpi_double_precision,0,
icomm,
ierrmpi)
629 call mpi_bcast(wy,ntheta_point,mpi_double_precision,0,
icomm,
ierrmpi)
630 call mpi_bcast(aphi,nphi_point,mpi_double_precision,0,
icomm,
ierrmpi)
631 call mpi_bcast(wphi,nphi_point,mpi_double_precision,0,
icomm,
ierrmpi)
644 real(8),
intent(in) :: xlow, xhi
645 integer,
intent(in) :: n
646 real(8),
intent(out) :: x(n), w(n)
649 real(8) :: p1, p2, p3, pp, xl, xm, z, z1
650 real(8),
parameter :: error=3.0
d-14
654 xm = 0.5d0*(xhi + xlow)
655 xl = 0.5d0*(xhi - xlow)
658 z = cos( dpi * (i - 0.25d0)/(n + 0.5d0) )
661 do while (abs(z1 - z) > error)
668 p1 = ( (2.0d0*j - 1.0d0)*z*p2 - (j - 1.0d0)*p3 )/j
671 pp = n*(z*p1 - p2) / (z*z - 1.0d0)
678 w(i) = 2.0d0*xl / ((1.0d0 - z*z) * pp*pp)
Module to include CAK radiation line force in (magneto)hydrodynamic models Computes both the force fr...
real(8), public gayley_qbar
real(8), public gayley_q0
logical cak_split
To treat source term in split or unsplit (default) fashion.
subroutine cak_init(phys_gamma)
Initialize the module.
subroutine cak_params_read(files)
Public method.
subroutine cak_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
Check time step for total radiation contribution.
subroutine gauss_legendre_quadrature(xlow, xhi, n, x, w)
Fast Gauss-Legendre N-point quadrature algorithm by G. Rybicki.
real(8), public cak_alpha
Line-ensemble parameters in the Gayley (1995) formalism.
subroutine cak_add_source(qdt, ixil, ixol, wct, w, x, energy, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
subroutine, public set_cak_force_norm(rstar, twind)
Compute some (unitless) variables for CAK force normalisation.
logical fix_vector_force_1d
To activate the pure radial vector CAK line force computation.
integer gcak1_
Extra slots to store quantities in w-array.
logical cak_vector_force
To activate the vector CAK line force computation.
subroutine get_velocity_gradient(ixil, ixol, vfield, x, idir, grad_vn)
Compute velocity gradient in direction 'idir' on a non-uniform grid.
integer cak_1d_opt
Switch to choose between the 1-D CAK line force options.
subroutine get_cak_force_radial(ixil, ixol, wct, w, x, gcak)
1-D CAK line force in the Gayley line-ensemble distribution parametrisation
subroutine get_gelectron(ixil, ixol, w, x, ge)
Compute continuum radiation force from Thomson scattering.
subroutine get_cak_force_vector(ixil, ixol, wct, w, x, gcak)
Vector CAK line force in the Gayley line-ensemble distribution parametrisation.
logical cak_1d_force
To activate the original CAK 1-D line force computation.
integer nthetaray
Amount of rays in radiation polar and radiation azimuthal direction.
subroutine rays_init(ntheta_point, nphi_point)
Initialise (theta',phi') radiation angles coming from stellar disc.
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter dpi
Pi.
double precision, parameter const_c
universal constants as specified in cgs units
double precision, parameter const_sigma
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_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
integer, parameter unitpar
file handle for IO
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision unit_length
Physical scaling factor for length.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision courantpar
The Courant (CFL) number used for the simulation.
integer ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
double precision unit_velocity
Physical scaling factor for velocity.
double precision unit_temperature
Physical scaling factor for temperature.
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_get_pthermal), pointer phys_get_pthermal
character(len=name_len) physics_type
String describing the physics type of the simulation.
Module with all the methods that users can customize in AMRVAC.