scal_schw Derived Type

type, public, extends(cpde_equation) :: scal_schw

A class that defines the evolution equations for a scalar field produced by a point charge in Schwarzschild coordinates decomposed into spherical harmonics.


Inherits

type~~scal_schw~~InheritsGraph type~scal_schw scal_schw type~cgf cgf type~scal_schw->type~cgf eq_flux_data type~rgf rgf type~scal_schw->type~rgf eq_coeffs, eq_lcoeffs, r_schw, r_star type~ref_element ref_element type~scal_schw->type~ref_element refelem type~rgfb rgfb type~scal_schw->type~rgfb eq_lambda, eq_s, eq_sinv type~scal_schw_eff scal_schw_eff type~scal_schw->type~scal_schw_eff effs type~cpde_equation cpde_equation type~scal_schw->type~cpde_equation type~tdc tdc type~scal_schw->type~tdc time_dep_coord type~element_cdata element_cdata type~cgf->type~element_cdata elems type~gf gf type~cgf->type~gf type~rgf->type~gf type~element_rdata element_rdata type~rgf->type~element_rdata elems type~element_boundary_rdata element_boundary_rdata type~rgfb->type~element_boundary_rdata elems type~gfb gfb type~rgfb->type~gfb type~eff_source eff_source type~scal_schw_eff->type~eff_source type~cpde_equation->type~cgf eq_data, eq_rhs_data, eq_tmp_data type~cgf_pointer cgf_pointer type~cpde_equation->type~cgf_pointer data_pointer type~equation equation type~cpde_equation->type~equation type~tdc->type~rgf dxdlambda, dxdxi, d2xdlambda2, d2xdxi2, d2xdlambdadxi, rm2m type~tdc->type~rgfb dxdlambda_b, dxdxi_b type~element_data element_data type~element_cdata->type~element_data type~cgf_pointer->type~cgf p type~gf_pointer gf_pointer type~cgf_pointer->type~gf_pointer type~element_boundary_data element_boundary_data type~element_boundary_rdata->type~element_boundary_data type~eff_source->type~cgf source type~element_rdata->type~element_data

Inherited by

type~~scal_schw~~InheritedByGraph type~scal_schw scal_schw type~sf_observer sf_observer type~sf_observer->type~scal_schw p

Contents

Source Code


Components

TypeVisibility AttributesNameInitial
integer(kind=ip), public :: ntmp

The number of temporary storage variables to allocate.

character(len=:), public, allocatable:: ename

The name of the system of equations

integer(kind=ip), public :: nmodes

The number of spherical harmonic modes. Should be moved to the specific implementation of the equation as this class should not be limited to systems using spherical harmonic decomposition.

integer(kind=ip), public :: nvars =3

The number of variables per mode. Should probably be changed to be the total number of variables in order to make this independent of spherical harmonic decomposition. And the value should be set by the equation implementation.

type(cgf), public, dimension(:,:), allocatable:: eq_data

A 2d array of complex grid functions that contains the data variables for the equation system. Should probably be changed to a 1d array to make this independent of spherical harmonic decomposition.

type(cgf), public, dimension(:,:), allocatable:: eq_rhs_data

A 2d array of complex grid functions that contains the rhs variables for the equation system. Should probably be changed to a 1d array to make this independent of spherical harmonic decomposition.

type(cgf), public, dimension(:,:,:), allocatable:: eq_tmp_data

A 3d array of complex grid functions that contains the temporary variables needed by the time integrator. Should probably be changed to a 2d array to make this independent of spherical harmonic decomposition.

type(cgf_pointer), public, dimension(:,:,:), allocatable:: data_pointer

A 3d array of pointers to complex grid functions that points to the RHS variables (:,:,-1), the data variables (:,:,0) and the temporary variables (:,:,1:ntmp) for the equation system. Should probably be changed to a 2d array to make this independent of spherical harmonic decomposition.

integer(kind=c_int), public, dimension(:), allocatable:: ll

A 1d array of c_int containing the l values of all the evolved modes.

integer(kind=c_int), public, dimension(:), allocatable:: mm

A 1d array of c_int containing the m values of all the evolved modes.

type(ref_element), public :: refelem

The ref_element contains the derivative matrix needed for approximating partial derivatives as well as the routines needed for calculating characteristic fluxes.

type(cgf), public, dimension(:,:), allocatable:: eq_flux_data

A 2d array of complex grid functions that stores flux data for the system.

type(rgf), public, dimension(:), allocatable:: eq_coeffs

A 1d array of real grid functions that stores the equation coefficients.

type(rgf), public, dimension(:), allocatable:: eq_lcoeffs

A 1d array of real grid functions that stores the l-dependent potential.

type(rgfb), public, dimension(:), allocatable:: eq_lambda

A 1d array of real boundary grid functions that stores the characteristic speeds at the element boundaries.

type(rgfb), public, dimension(:,:), allocatable:: eq_s

A 2d array of real boundary grid functions that stores the matrix used to convert from characteristic to evolved variables.

type(rgfb), public, dimension(:,:), allocatable:: eq_sinv

A 2d array of real boundary grid functions that stores the matrix used to convert from evolve to characteristic variables.

type(rgf), public :: r_schw

A real grid function that stores the Schwarzschild radial coordinate.

type(rgf), public :: r_star

A real grid function that stores the tortoise radial coordinate.

type(scal_schw_eff), public :: effs

The effective source for a scalar point charge on generic orbits in a Schwarzschild spacetime.

type(tdc), public :: time_dep_coord

A time dependent coordinate class for the scalar wave equation in Tortoise coordinates.


Type-Bound Procedures

procedure, public :: set_to_zero => cpde_set_to_zero

The set_to_zero routine is provided by cpde_set_to_zero.

  • interface

    public module subroutine cpde_set_to_zero(this, dest)

    The interface for the PDE version of set_to_zero. This interface is consistent with eq_set_to_zero_interface.

    Arguments

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

    The routine is called on this equation.

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

    Can be either -1 (RHS), 0 (main) or between 1 and ntmp (temporary).

procedure, public :: update_vars => cpde_update_vars

The update_vars routine is provided by cpde_update_vars.

  • interface

    public module subroutine cpde_update_vars(this, source, dest, source2, scalar, scalar2)

    The interface for the PDE version of update_vars. This interface is consistent with eq_update_vars_interface.

    Arguments

    Type IntentOptional AttributesName
    class(cpde_equation), intent(inout), target:: this

    The routine is called on this equation.

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

    The first source for the update. Can be either -1 (RHS), 0 (main) or between 1 and ntmp (temporary).

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

    The destination for the update. Can be either -1 (RHS), 0 (main) or between 1 and ntmp (temporary).

    integer(kind=ip), intent(in), optional :: source2

    The second source for the update. Can be either -1 (RHS), 0 (main) or between 1 and ntmp (temporary).

    real(kind=wp), intent(in), optional :: scalar

    The scalar multiplying the first source.

    real(kind=wp), intent(in), optional :: scalar2

    The scalar multiplying the second source.

procedure, public :: output => cpde_output

The output routine is provided by cpde_output.

  • interface

    public module subroutine cpde_output(this)

    The interface for the PDE version of output. This interface is consistent with eq_output.

    Arguments

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

    The routine is called on this equation.

procedure, public :: print_data => cpde_print_data

The print_data routine is provided by cpde_print_data.

  • interface

    public module subroutine cpde_print_data(this)

    The interface for the PDE version of print_data. This interface is consistent with eq_print_data.

    Arguments

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

    The routine is called on this equation.

procedure, public :: init => scal_schw_init

The init routine is provided by scal_schw_init.

  • interface

    public module subroutine scal_schw_init(this)

    The interface for the scal_schw version of init. This interface is consistent with eq_init_interface.

    Arguments

    Type IntentOptional AttributesName
    class(scal_schw), intent(inout), target:: this

    The equation that is being initialized.

procedure, public :: rhs => scal_schw_rhs

The rhs routine is provided by scal_schw_rhs.

  • interface

    public module subroutine scal_schw_rhs(this)

    The interface for the scal_schw version of rhs. This interface is consistent with eq_rhs_interface.

    Arguments

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

    The equation for which the RHS is calculated.

procedure, public :: save_globals_1 => scal_schw_save_globals_1

The save_globals_1 routine is provided by scal_schw_save_globals_1.

procedure, public :: save_globals_2 => scal_schw_save_globals_2

The save_globals_2 routine is provided by scal_schw_save_globals_2.

procedure, public :: load_globals => scal_schw_load_globals

The load_globals routine is provided by scal_schw_load_globals.

  • interface

    public module subroutine scal_schw_load_globals(this)

    The interface for the scal_schw version of load_globals. This interface is consistent with eq_load_globals.

    Arguments

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

    The routine is called on this equation.

generic, public :: flux => scal_schw_flux

A routine to calculate fluxes is provided by scal_schw_flux under the generic name flux.

procedure, public :: scal_schw_flux

Provides a routine to calculate fluxes.

  • interface

    public module subroutine scal_schw_flux(this)

    The interface for the scal_schw_flux routine for calculating fluxes. Can be called under the generic name flux.

    Arguments

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

    The routine is called on this equation.

procedure, public :: read_all_modes

Provides a routine to read in external initial data.

  • interface

    public module subroutine read_all_modes(this)

    The interface for the read_all_modes routine that reads external initial data.

    Arguments

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

    The routine is called on this equation.

procedure, public :: output_coords

Provides a routine to output the coordinates in a format needed for a frequency domain python code that can provide initial data.

  • interface

    public module subroutine output_coords(this)

    The interface for the output_coords routine that writes out the coordinates for an external python frequency domain code that calculates initial data.

    Arguments

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

    The routine is called on this equation.

procedure, public :: tortoise_to_hyperboloidal

Provides a routine to transform from tortoise to hyperboloid that is used when reading in external data.

  • interface

    public module subroutine tortoise_to_hyperboloidal(this, elem, node, inner, dpsidt, dpsidr)

    The interface for the tortoise_to_hyperboloidal routine that converts from Tortoise to hyperboloidal coordinates when reading in external initial data.

    Arguments

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

    The routine is called on this equation.

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

    The index of the element where the transformation is calculated.

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

    The index of the node within that element where the transformation is calculated.

    logical, intent(in) :: inner

    If .true. use the transformation for the inner region. If .false. use the transformation for the outer region.

    complex(kind=wp), intent(inout) :: dpsidt

    On input holds the time derivative in Tortoise coordinates. On output holds the time derivative in hyperboloidal coordinates.

    complex(kind=wp), intent(inout) :: dpsidr

    On input holds the radial derivative in Tortoise coordinates. On output holds the radial derivative in hyperboloidal coordinates.

Source Code

  type, extends(cpde_equation) :: scal_schw
  !! A class that defines the evolution equations for a scalar field produced
  !! by a point charge in Schwarzschild coordinates decomposed into
  !! spherical harmonics.
    integer(c_int), dimension(:), allocatable :: ll
    !! A 1d array of c_int containing the l values of all the evolved modes.
    integer(c_int), dimension(:), allocatable :: mm
    !! A 1d array of c_int containing the m values of all the evolved modes.
    type(ref_element) :: refelem
    !! The [[ref_element]] contains the derivative matrix needed for
    !! approximating partial derivatives as well as the routines needed
    !! for calculating characteristic fluxes.
    type(cgf), dimension(:,:), allocatable :: eq_flux_data
    !! A 2d array of complex grid functions that stores flux data for the
    !! system.
    type(rgf), dimension(:), allocatable :: eq_coeffs
    !! A 1d array of real grid functions that stores the equation coefficients.
    type(rgf), dimension(:), allocatable :: eq_lcoeffs
    !! A 1d array of real grid functions that stores the l-dependent
    !! potential.
    type(rgfb), dimension(:), allocatable :: eq_lambda
    !! A 1d array of real boundary grid functions that stores the
    !! characteristic speeds at the element boundaries.
    type(rgfb), dimension(:,:), allocatable :: eq_s
    !! A 2d array of real boundary grid functions that stores the
    !! matrix used to convert from characteristic to evolved variables.
    type(rgfb), dimension(:,:), allocatable :: eq_sinv
    !! A 2d array of real boundary grid functions that stores the
    !! matrix used to convert from evolve to characteristic variables.
    type(rgf) :: r_schw
    !! A real grid function that stores the Schwarzschild radial coordinate.
    type(rgf) :: r_star
    !! A real grid function that stores the tortoise radial coordinate.
    type(scal_schw_eff) :: effs
    !! The effective source for a scalar point charge on generic orbits in a
    !! Schwarzschild spacetime.
    type(tdc) :: time_dep_coord
    !! A time dependent coordinate class for the scalar wave equation in
    !! Tortoise coordinates.
  contains
    procedure :: init => scal_schw_init
    !! The [[equation:init]] routine is provided by [[scal_schw_init]].
    procedure :: rhs => scal_schw_rhs
    !! The [[equation:rhs]] routine is provided by [[scal_schw_rhs]].
    procedure :: save_globals_1 => scal_schw_save_globals_1
    !! The [[equation:save_globals_1]] routine is provided by
    !! [[scal_schw_save_globals_1]].
    procedure :: save_globals_2 => scal_schw_save_globals_2
    !! The [[equation:save_globals_2]] routine is provided by
    !! [[scal_schw_save_globals_2]].
    procedure :: load_globals => scal_schw_load_globals
    !! The [[equation:load_globals]] routine is provided by
    !! [[scal_schw_load_globals]].
    generic :: flux => scal_schw_flux
    !! A routine to calculate fluxes is provided by [[scal_schw_flux]] under
    !! the generic name [[scal_schw:flux]].
    procedure :: scal_schw_flux
    !! Provides a routine to calculate fluxes.
    procedure :: read_all_modes
    !! Provides a routine to read in external initial data.
    procedure :: output_coords
    !! Provides a routine to output the coordinates in a format needed for
    !! a frequency domain python code that can provide initial data.
    procedure :: tortoise_to_hyperboloidal
    !! Provides a routine to transform from tortoise to hyperboloid that
    !! is used when reading in external data.
  end type scal_schw