scalar_schw_eff_implementation Submodule

The implementation of the interfaces provided in scalar_schw_eff.


Uses

  • module~~scalar_schw_eff_implementation~~UsesGraph module~scalar_schw_eff_implementation scalar_schw_eff_implementation module~scalar_schw_eff scalar_schw_eff module~scalar_schw_eff_implementation->module~scalar_schw_eff module~kinds kinds module~scalar_schw_eff->module~kinds module~effective_source effective_source module~scalar_schw_eff->module~effective_source iso_c_binding iso_c_binding module~scalar_schw_eff->iso_c_binding module~effective_source->module~kinds module~effective_source->iso_c_binding module~grid_function grid_function module~effective_source->module~grid_function module~world_tube world_tube module~effective_source->module~world_tube module~grid_function->module~kinds module~element element module~grid_function->module~element module~world_tube->module~kinds module~world_tube->module~grid_function module~element->module~kinds

Contents


Interfaces

interface

  • subroutine initialize_source(nmodes, l, m, mass) bind(c, name="0")

    Interface to the C++ routine init_source for initializing the C++ effective source class.

    Arguments

    Type IntentOptional AttributesName
    integer(kind=c_int) :: nmodes

    The number of modes.

    integer(kind=c_int), dimension(nmodes):: l

    A 1d array that contains the -values.

    integer(kind=c_int), dimension(nmodes):: m

    A 1d array that contains the -values.

    real(kind=c_double) :: mass

    The mass of the black hole.

interface

  • subroutine set_time_window_coeffs(tfac, dtfac_dt, d2tfac_dt2, do_smooth_after_lmax, nmodes) bind(c, name="0")

    Interface to the C++ routine set_time_window that sets the time window information in the C++ effective source class.

    Arguments

    Type IntentOptional AttributesName
    real(kind=c_double) :: tfac

    The time window factor.

    real(kind=c_double) :: dtfac_dt

    The first time derivative of the time window factor.

    real(kind=c_double) :: d2tfac_dt2

    The second time derivative of the time window factor.

    integer(kind=c_int) :: do_smooth_after_lmax

    Only set the time window for do_smooth_after_lmax.

    integer(kind=c_int) :: nmodes

    The number of modes.

interface

  • subroutine set_particle_pos(r, phi, ur, En, Lz, ar, aphi, dardt, daphidt, d2ardt2, d2aphidt2, nmodes) bind(c, name="0")

    Interface to the C++ routine set_particle that sets the particle information in the C++ effective source class.

    Arguments

    Type IntentOptional AttributesName
    real(kind=c_double) :: r

    The radial coordinate of the particle (in Schwarzschild coordinates), .

    real(kind=c_double) :: phi

    The azimuthal angle of the particle, .

    real(kind=c_double) :: ur

    The radial component of the 4-velocity, .

    real(kind=c_double) :: En

    The energy per unit mass of the particle, .

    real(kind=c_double) :: Lz

    The angular momentum per unit mass of the particle .

    real(kind=c_double) :: ar

    The radial component of the 4-acceleration, .

    real(kind=c_double) :: aphi

    The -component of the 4-acceleration, .

    real(kind=c_double) :: dardt

    The derivative of with respect to Schwarzschild coordinate time.

    real(kind=c_double) :: daphidt

    The derivative of with respect to Schwarzschild coordinate time.

    real(kind=c_double) :: d2ardt2

    The second derivative of with respect to Schwarzschild coordinate time.

    real(kind=c_double) :: d2aphidt2

    The second derivative of with respect to Schwarzschild coordinate time.

    integer(kind=c_int) :: nmodes

    The number of modes.

interface

  • subroutine evaluate_source_all(mode, n, r, win, dwin, d2win, sre, sim) bind(c, name="0")

    Interface to the C++ routine eval_source_all that evaluates the effective source on the locations provided in a 1d input array.

    Arguments

    Type IntentOptional AttributesName
    integer(kind=c_int) :: mode

    The mode.

    integer(kind=c_int) :: n

    The number of elements in the input radial coordinate array, .

    real(kind=c_double) :: r(*)

    The radial coordinates in Schwarzschild coordinates, .

    real(kind=c_double) :: win(*)

    The window function, .

    real(kind=c_double) :: dwin(*)

    The radial derivative of the window function, .

    real(kind=c_double) :: d2win(*)

    The second radial derivative of the window function, .

    real(kind=c_double) :: sre(*)

    The real part of the effective source, .

    real(kind=c_double) :: sim(*)

    The imaginary part of the effective source, .

interface

  • subroutine get_phi(mode, r, phire, phiim) bind(c, name="0")

    Interface to the C++ routine Phi that evaluates the singuler field at a given radial coordinate in Schwarzschild coordinates.

    Arguments

    Type IntentOptional AttributesName
    integer(kind=c_int) :: mode

    The mode.

    real(kind=c_double) :: r

    The radial coordinate, , in Schwarzschild coordinates.

    real(kind=c_double) :: phire

    On output the real part of the singular field .

    real(kind=c_double) :: phiim

    On output the imaginary part of the singular field .

interface

  • subroutine get_dphidt(mode, r, dphidtre, dphidtim) bind(c, name="0")

    Interface to the C++ routine dPhi_dt that evaluates the time derivative of the singuler field at a given radial coordinate in Schwarzschild coordinates.

    Arguments

    Type IntentOptional AttributesName
    integer(kind=c_int) :: mode

    The mode.

    real(kind=c_double) :: r

    The radial coordinate, , in Schwarzschild coordinates.

    real(kind=c_double) :: dphidtre

    On output the real part of the time derivative of the singular field .

    real(kind=c_double) :: dphidtim

    On output the imaginary part of the time derivative of the singular field .

interface

  • subroutine get_dphidr(mode, r, dphidrre, dphidrim) bind(c, name="0")

    Interface to the C++ routine dPhi_dt that evaluates the radial derivative of the singuler field at a given radial coordinate in Schwarzschild coordinates.

    Arguments

    Type IntentOptional AttributesName
    integer(kind=c_int) :: mode

    The mode.

    real(kind=c_double) :: r

    The radial coordinate, , in Schwarzschild coordinates.

    real(kind=c_double) :: dphidrre

    On output the real part of the radial derivative of the singular field .

    real(kind=c_double) :: dphidrim

    On output the imaginary part of the radial derivative of the singular field .


Module Procedures

module procedure scal_schw_eff_init module subroutine scal_schw_eff_init(this, nmodes, nvars, l, m, mass)

Arguments

Type IntentOptional AttributesName
class(scal_schw_eff), intent(inout) :: this

On return, the constructed object.

integer(kind=ip), intent(in) :: nmodes

The number of modes.

integer(kind=ip), intent(in) :: nvars

The number of variables.

integer(kind=c_int), intent(in), dimension(nmodes):: l

A 1d array of size nmodes containing the -values of the modes.

integer(kind=c_int), intent(in), dimension(nmodes):: m

A 1d array of size nmodes containing the m-values of the modes.

real(kind=wp), intent(in) :: mass

The mass of the black hole.

module procedure scal_schw_eff_set_particle_pos module subroutine scal_schw_eff_set_particle_pos(this, r, phi, ur, En, Lz, ar, aphi, dardt, daphidt, d2ardt2, d2aphidt2)

Arguments

Type IntentOptional AttributesName
class(scal_schw_eff), intent(inout) :: this

The routine is called on this object.

real(kind=wp), intent(in) :: r

The radial coordinate of the particle (in Schwarzschild coordinates), .

real(kind=wp), intent(in) :: phi

The azimuthal angle of the particle, .

real(kind=wp), intent(in) :: ur

The radial component of the 4-velocity, .

real(kind=wp), intent(in) :: En

The energy per unit mass of the particle, .

real(kind=wp), intent(in) :: Lz

The angular momentum per unit mass of the particle .

real(kind=wp), intent(in) :: ar

The radial component of the 4-acceleration, .

real(kind=wp), intent(in) :: aphi

The -component of the 4-acceleration, .

real(kind=wp), intent(in) :: dardt

The derivative of with respect to Schwarzschild coordinate time.

real(kind=wp), intent(in) :: daphidt

The derivative of with respect to Schwarzschild coordinate time.

real(kind=wp), intent(in) :: d2ardt2

The second derivative of with respect to Schwarzschild coordinate time.

real(kind=wp), intent(in) :: d2aphidt2

The second derivative of with respect to Schwarzschild coordinate time.

module procedure scal_schw_eff_set_time_window module subroutine scal_schw_eff_set_time_window(this, tfac, dtfac_dt, d2tfac_dt2, do_smooth_after_lmax)

Arguments

Type IntentOptional AttributesName
class(scal_schw_eff), intent(inout) :: this

The routine is called on this object.

real(kind=wp), intent(in) :: tfac

The current value for the time window.

real(kind=wp), intent(in) :: dtfac_dt

The current value for the time derivative of the time window.

real(kind=wp), intent(in) :: d2tfac_dt2

The current value for the second time derivative of the time window.

integer(kind=ip), intent(in) :: do_smooth_after_lmax

If do_smooth_after_lmax use tfac=1, dtfac_dt=0, d2tfac_dt2=0. This allows for smooth turn on of the effective source for the modes with do_smooth_after_lmax, while the rest gets turned on instantaneously (e.g. when external initial data is available).

module procedure scal_schw_eff_evaluate_source module subroutine scal_schw_eff_evaluate_source(this, r, wt)

Arguments

Type IntentOptional AttributesName
class(scal_schw_eff), intent(inout) :: this

The routine is called on this object.

type(rgf), intent(in) :: r

A real values grid function that contain the radial coordinate (in Schwarzschild coordinates).

type(wtube), intent(in) :: wt

A world-tube object that ensures that the effective source is only non-zero inside the world-tube.

module procedure scal_schw_eff_get_singular module subroutine scal_schw_eff_get_singular(this, r, mode, psi)

Arguments

Type IntentOptional AttributesName
class(scal_schw_eff), intent(inout) :: this

The routine is called on this object.

real(kind=wp), intent(in) :: r

The radial coordinate (Schwarzschild coordinates).

integer(kind=ip), intent(in) :: mode

The mode.

complex(kind=wp), intent(out), dimension(:):: psi

A 1d-array of size nvars of complex values that on return contains the singular field for all variables.

module procedure scal_schw_eff_get_dsingular_dt module subroutine scal_schw_eff_get_dsingular_dt(this, r, mode, dpsidt)

Arguments

Type IntentOptional AttributesName
class(scal_schw_eff), intent(inout) :: this

The routine is called on this object.

real(kind=wp), intent(in) :: r

The radial coordinate (Schwarzschild coordinates).

integer(kind=ip), intent(in) :: mode

The mode.

complex(kind=wp), intent(out), dimension(:):: dpsidt

A 1d-array of size nvars of complex values that on return contains the time derivative of the singular field for all variables.

module procedure scal_schw_eff_get_dsingular_dr module subroutine scal_schw_eff_get_dsingular_dr(this, r, mode, dpsidr)

Arguments

Type IntentOptional AttributesName
class(scal_schw_eff), intent(inout) :: this

The routine is called on this object.

real(kind=wp), intent(in) :: r

The radial coordinate (Schwarzschild coordinates).

integer(kind=ip), intent(in) :: mode

The mode.

complex(kind=wp), intent(out), dimension(:):: dpsidr

A 1d-array of size nvars of complex values that on return contains the radial derivative of the singular field for all variables.