bssn Derived Type

type, public, extends(tpo) :: bssn

TYPE representing the \(\mathrm{BSSNOK}\) formulation of the Einstein equations


Inherits

type~~bssn~~InheritsGraph type~bssn bssn grid_function grid_function type~bssn->grid_function Gamma_u, A_BSSN3_ll, g_BSSN3_ll, Ricci_ll, GC, GC_parts grid_function_scalar grid_function_scalar type~bssn->grid_function_scalar phi, trK, Ricci_scalar timer timer type~bssn->timer bssn_computer_timer type~tpo tpo type~bssn->type~tpo type~tpo->grid_function coords, shift_u, g_phys3_ll, K_phys3_ll, S, MC, MC_parts, S_parts type~tpo->grid_function_scalar rad_coord, lapse, HC, HC_parts, rho, rho_parts type~tpo->timer grid_timer, importer_timer level level type~tpo->level levels

Contents

Source Code


Components

Type Visibility Attributes Name Initial
type(grid_function), public :: A_BSSN3_ll

Conformal traceless extrinsic curvature

type(grid_function), public :: GC

Connection constraint computed with the \(\mathrm{ID}\) on the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: GC_int

Integral of the connection constraint computed with the \(\mathrm{ID}\) on the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: GC_l2

norm of the connection constraint computed with the \(\mathrm{ID}\) on the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: GC_loo

norm, i.e., supremum of the absolute value, of the connection constraint computed with the \(\mathrm{ID}\) on the mesh

type(grid_function), public :: GC_parts

Connection constraint computed with the \(\mathrm{BSSNOK}\) \(\mathrm{ID}\) on the mesh, and the hydrodynamical \(\mathrm{ID}\) mapped from the particles to the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: GC_parts_int

Integral of the connection constraint computed with the \(\mathrm{BSSNOK}\) \(\mathrm{ID}\) on the mesh, and the hydrodynamical \(\mathrm{ID}\) mapped from the particles to the mesh with the \(\mathrm{ID}\) on the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: GC_parts_l2

norm of the connection constraint computed with the \(\mathrm{BSSNOK}\) \(\mathrm{ID}\) on the mesh, and the hydrodynamical \(\mathrm{ID}\) mapped from the particles to the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: GC_parts_loo

norm, i.e., supremum of the absolute value, of the connection constraint computed with the \(\mathrm{BSSNOK}\) \(\mathrm{ID}\) on the mesh, and the hydrodynamical \(\mathrm{ID}\) mapped from the particles to the mesh

type(grid_function), public :: Gamma_u

Conformal connection

type(grid_function_scalar), public :: HC

Grid scalar function storing the Hamiltonian constraint (violations) computed using the ID on the mesh

double precision, public, DIMENSION(:), ALLOCATABLE :: HC_int

Integral the Hamiltonian constraint computed using the stress-energy tensor mapped from the particles

double precision, public, DIMENSION(:), ALLOCATABLE :: HC_l2

norm of the Hamiltonian constraint computed on the mesh

double precision, public, DIMENSION(:), ALLOCATABLE :: HC_loo

norm of the Hamiltonian constraint computed on the mesh (i.e., its maximum)

type(grid_function_scalar), public :: HC_parts

Grid scalar function storing the Hamiltonian constraint (violations) computed using the stress-energy tensor mapped from the particles to the mesh

double precision, public, DIMENSION(:), ALLOCATABLE :: HC_parts_int

Integral the Hamiltonian constraint computed on the mesh

double precision, public, DIMENSION(:), ALLOCATABLE :: HC_parts_l2

norm of the Hamiltonian constraint computed on the mesh, using the stress-energy tensor mapped from the particles

double precision, public, DIMENSION(:), ALLOCATABLE :: HC_parts_loo

norm of the Hamiltonian constraint computed on the mesh, using the stress-energy tensor mapped from the particles (i.e., its maximum)

type(grid_function), public :: K_phys3_ll

Grid function storing the extrinsic curvature

type(grid_function), public :: MC

Grid function storing the momentum constraint (violations) computed using the ID on the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: MC_int

Integral of the momentum constraint computed on the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: MC_l2

norm of the momentum constraint computed on the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: MC_loo

norm of the momentum constraint computed on the mesh (i.e., its maximum)

type(grid_function), public :: MC_parts

Grid function storing the momentum constraint (violations) computed using the stress-energy tensor mapped from the particles to the mesh

double precision, public, DIMENSION(:,:), ALLOCATABLE :: MC_parts_int

Integral of the momentum constraint computed using the stress-energy tensor mapped from the particles

double precision, public, DIMENSION(:,:), ALLOCATABLE :: MC_parts_l2

norm of the Hamiltonian constraint computed on the mesh, using the stress-energy tensor mapped from the particles

double precision, public, DIMENSION(:,:), ALLOCATABLE :: MC_parts_loo

norm of the momentum constraint computed on the mesh, using the stress-energy tensor mapped from the particles (i.e., its maximum)

type(grid_function), public :: Ricci_ll

Ricci tensor

type(grid_function_scalar), public :: Ricci_scalar

Ricci scalar

type(grid_function), public :: S

Grid function storing the matter source in the momentum constraint computed using the ID on the mesh, multiplied by

type(grid_function), public :: S_parts

Grid function storing the matter source in the momentum constraint computed using the stress-energy tensor mapped from the particles to the mesh, multiplied by

type(timer), public :: bssn_computer_timer

Timer that times how long it takes to compute the \(\mathrm{BSSNOK}\) variables on the refined mesh

integer, public :: call_flag = 0

Flag set to a value different than 0 if the SUBROUTINE compute_and_print_bssn_variables is called

integer, public :: cons_step

Constraint violations are printed to file every cons_step-th grid point

type(grid_function), public :: coords

Grid function storing the Cartesian coordinates

logical, public :: export_bin

.TRUE. if the binary files for SPHINCS_BSSN are to be exported, .FALSE. otherwise

logical, public :: export_constraints

.TRUE. if the constraint violations are to be printed to file, .FALSE. otherwise

logical, public :: export_constraints_details

.TRUE. if the points at which the constraints violations are within the intervals , with and , are to be printed to file; .FALSE. otherwise

logical, public :: export_constraints_x

.TRUE. if only the constrain violations on the axis are to be

logical, public :: export_constraints_xy

.TRUE. if only the constrain violations on the plane are to be

logical, public :: export_form_x

.TRUE. if the \(\mathrm{ID}\) in the formatted files is to be on the x axis only, .FALSE. otherwise

logical, public :: export_form_xy

.TRUE. if the \(\mathrm{ID}\) in the formatted files is to be on the xy plane only, .FALSE. otherwise

type(grid_function), public :: g_BSSN3_ll

Conformal spatial metric

type(grid_function), public :: g_phys3_ll

Grid function storing the spatial metric

type(timer), public :: grid_timer

Timer that times how long it takes to set up the grid and allocate the grid functions

type(timer), public :: importer_timer

Timer that times how long it takes to import the \(\texttt{LORENE}\) ID on the mesh

type(grid_function_scalar), public :: lapse

Grid scalar function storing the lapse function

type(level), public, DIMENSION(:), ALLOCATABLE :: levels

Array containing the information on each refinement level

integer, public :: n_matter

Number of matter objects in the physical system

integer, public :: nlevels

Number of refinement levels

integer, public, DIMENSION(:), ALLOCATABLE :: npoints_xaxis

Array containing the number of mesh points of the highest-resolution refinement level across the x-axis-size of the matter objects

type(grid_function_scalar), public :: phi

Conformal factor

type(grid_function_scalar), public :: rad_coord

Grid scalar function storing the radial coordinates of each grid point

type(grid_function_scalar), public :: rho

Grid scalar function storing the matter source in the Hamiltonian constraint computed using the ID on the mesh, multiplied by

type(grid_function_scalar), public :: rho_parts

Grid scalar function storing the matter source in the Hamiltonian constraint computed using the stress-energy tensor mapped from the particles to the mesh, multiplied by

type(grid_function), public :: shift_u

Grid function storing the shift vector

integer, public :: tpo_id_number

Negative integer that identifies the tpo object

type(grid_function_scalar), public :: trK

Trace of extrinsic curvature


Constructor

public interface bssn


Finalization Procedures

final :: destructor

Destructor; finalizes members from both TYPES tpo and bssn, by calling destruct_bssn


Type-Bound Procedures

procedure, public :: abs_values_in

  • interface

    public module subroutine abs_values_in(this, lower_bound, upper_bound, constraint, l, export, unit_analysis, cnt)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    double precision, intent(in) :: lower_bound
    double precision, intent(in) :: upper_bound
    double precision, intent(in), DIMENSION(:,:,:) :: constraint
    integer, intent(in) :: l
    logical, intent(in) :: export
    integer, intent(in) :: unit_analysis
    integer, intent(out) :: cnt

procedure, public, NON_OVERRIDABLE :: analyze_constraint

Analyze a constraint (or an arbitrary scalar grid function) by examining its values at the refined mesh. Computes the number of mesh points at which the scalar grid function has values lying within predefined (hard-coded) intervals

  • interface

    public module subroutine analyze_constraint(this, l, constraint, name_constraint, unit_logfile, name_analysis, l2_norm, loo_norm, integral, source)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l
    double precision, intent(in), DIMENSION(:,:,:) :: constraint
    character(len=*), intent(in) :: name_constraint
    integer, intent(in) :: unit_logfile
    character(len=*), intent(in) :: name_analysis
    double precision, intent(out) :: l2_norm
    double precision, intent(out) :: loo_norm
    double precision, intent(out) :: integral
    double precision, intent(in), optional, DIMENSION(:,:,:) :: source

procedure, public, NON_OVERRIDABLE :: compute_adm_momentum_fluid_m2p

Computes an estimate of the linear momentum of the fluid using the \(\mathrm{SPH}\) hydro fields, and the spacetime metric mapped from the mesh

  • interface

    public module subroutine compute_adm_momentum_fluid_m2p(this, parts, adm_mom)

    Computes an estimate of the linear momentum of the fluid using the \(\mathrm{SPH}\) hydro fields, and the spacetime metric mapped from the mesh

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(in) :: this

    tpo object which this PROCEDURE is a member of

    class(particles), intent(in) :: parts

    particles object used to map the metric from the mesh to the particles, and to call the recovery procedures

    double precision, intent(out), DIMENSION(3) :: adm_mom

    ADM linear momentum of the fluid computed using the metric mapped with the mesh-to-particle mapping

generic, public :: compute_and_print_tpo_constraints => compute_and_print_tpo_constraints_grid, compute_and_print_tpo_constraints_particles

Overloaded PROCEDURE to compute the constraints using only the \(\mathrm{ID}\) on the refined mesh, or the spacetime \(\mathrm{ID}\) on the refined mesh and the hydrodynamical \(\mathrm{ID}\) mapped from the particles to the refined mesh

procedure, public :: compute_and_print_tpo_constraints_grid => compute_and_print_bssn_constraints_grid

Computes the \(\mathrm{BSSNOK}\) constraints using the full \(\mathrm{ID}\) on the refined mesh, prints a summary with the statistics for the constraints. Optionally, prints the constraints to a formatted file to be read by , and analyze the constraints by calling analyze_constraint

  • interface

    public module subroutine compute_and_print_bssn_constraints_grid(this, id, namefile, name_logfile, points)

    Interface to compute_and_print_tpo_constraints_grid

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout), TARGET :: this

    bssn object to which this PROCEDURE is bound

    class(idbase), intent(inout) :: id

    idbase object used to read the hydrodynamical \(\mathrm{ID}\) to the mesh

    character(len=*), intent(inout) :: namefile
    character(len=*), intent(inout) :: name_logfile
    double precision, intent(in), optional, DIMENSION(:,:,:,:), TARGET :: points

procedure, public :: compute_and_print_tpo_constraints_particles => compute_and_print_bssn_constraints_particles

Computes the \(\mathrm{BSSNOK}\) constraints using the \(\mathrm{BSSNOK}\) \(\mathrm{ID}\) on the refined mesh and the hydrodynamical \(\mathrm{ID}\) mapped from the particles to the mesh, prints a summary with the statistics for the constraints. Optionally, prints the constraints to a formatted file to be read by , and analyze the constraints by calling analyze_constraint

  • interface

    public module subroutine compute_and_print_bssn_constraints_particles(this, parts_obj, namefile, name_logfile, points)

    Interface to compute_and_print_tpo_constraints_particles

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout), TARGET :: this

    bssn object to which this PROCEDURE is bound

    class(particles), intent(inout) :: parts_obj

    particles object used to map the hydrodynamical \(\mathrm{ID}\) to the mesh

    character(len=*), intent(inout) :: namefile
    character(len=*), intent(inout) :: name_logfile
    double precision, intent(in), optional, DIMENSION(:,:,:,:), TARGET :: points

procedure, public :: compute_and_print_tpo_variables => compute_and_print_bssn_variables

Computes the \(\mathrm{BSSNOK}\) variables at the particle positions, and optionally prints them to a binary file to be read by and , and to a formatted file to be read by , by calling print_formatted_id_tpo_variables

  • interface

    public module subroutine compute_and_print_bssn_variables(this, namefile)

    Interface to compute_and_print_tpo_variables

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout) :: this

    bssn object to which this PROCEDURE is bound

    character(len=*), intent(inout), optional :: namefile

procedure, public :: compute_ricci

Computes the Ricci tensor and the Ricci scalar on the mesh

  • interface

    public module subroutine compute_ricci(this)

    Computes the Ricci tensor and the Ricci scalar on the mesh

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout) :: this

    bssn object to which this PROCEDURE is bound

procedure, public :: deallocate_fields => deallocate_bssn_fields

Deallocates memory for the bssn member arrays

  • interface

    public module subroutine deallocate_bssn_fields(this)

    Interface to deallocate_fields

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout) :: this

    bssn object to which this PROCEDURE is bound

procedure, public, NON_OVERRIDABLE :: deallocate_standard_tpo_variables

Deallocates memory for the standard 3+1 fields

  • interface

    public module subroutine deallocate_standard_tpo_variables(this)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this

procedure, public :: define_allocate_fields => allocate_bssn_fields

Allocates memory for the bssn member grid functions

  • interface

    public module subroutine allocate_bssn_fields(this)

    Interface to define_allocate_fields

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout) :: this

    bssn object to which this PROCEDURE is bound

procedure, public :: destruct_bssn

Destructor for the EXTENDED TYPE bssn, not ABSTRACT TYPE tpo

  • interface

    public module subroutine destruct_bssn(this)

    Interface to destruct_bssn

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout) :: this

    bssn object to which this PROCEDURE is bound

procedure, public :: get_HC

  • interface

    public module function get_HC(this, i, j, k, l) result(HC_value)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: i
    integer, intent(in) :: j
    integer, intent(in) :: k
    integer, intent(in) :: l

    Return Value double precision

procedure, public :: get_HC_parts

  • interface

    public module function get_HC_parts(this, i, j, k, l) result(HC_value)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: i
    integer, intent(in) :: j
    integer, intent(in) :: k
    integer, intent(in) :: l

    Return Value double precision

procedure, public :: get_MC

  • interface

    public module function get_MC(this, i, j, k, l) result(MC_value)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: i
    integer, intent(in) :: j
    integer, intent(in) :: k
    integer, intent(in) :: l

    Return Value double precision, DIMENSION(3)

procedure, public :: get_MC_parts

  • interface

    public module function get_MC_parts(this, i, j, k, l) result(MC_value)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: i
    integer, intent(in) :: j
    integer, intent(in) :: k
    integer, intent(in) :: l

    Return Value double precision, DIMENSION(3)

procedure, public :: get_dx

  • interface

    public module function get_dx(this, l) result(dx)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value double precision

procedure, public :: get_dy

  • interface

    public module function get_dy(this, l) result(dy)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value double precision

procedure, public :: get_dz

  • interface

    public module function get_dz(this, l) result(dz)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value double precision

procedure, public :: get_grid_point

  • interface

    public module function get_grid_point(this, i, j, k, l) result(grid_point)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: i
    integer, intent(in) :: j
    integer, intent(in) :: k
    integer, intent(in) :: l

    Return Value double precision, DIMENSION(3)

procedure, public :: get_levels

  • interface

    public module function get_levels(this) result(levels)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this

    Return Value type(level), DIMENSION(:), ALLOCATABLE

procedure, public :: get_ngrid_x

  • interface

    public module function get_ngrid_x(this, l) result(ngrid_x)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value integer

procedure, public :: get_ngrid_y

  • interface

    public module function get_ngrid_y(this, l) result(ngrid_y)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value integer

procedure, public :: get_ngrid_z

  • interface

    public module function get_ngrid_z(this, l) result(ngrid_z)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value integer

procedure, public :: get_nlevels

  • interface

    public module function get_nlevels(this) result(nlevels)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this

    Return Value integer

procedure, public :: get_xR

  • interface

    public module function get_xR(this, l) result(xR)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value double precision

procedure, public :: get_yR

  • interface

    public module function get_yR(this, l) result(yR)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value double precision

procedure, public :: get_zR

  • interface

    public module function get_zR(this, l) result(zR)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    integer, intent(in) :: l

    Return Value double precision

procedure, public :: print_formatted_id_tpo_variables => print_formatted_id_bssn_variables

Prints the \(\mathrm{BSSNOK}\) \(\mathrm{ID}\) to a formatted file

  • interface

    public module subroutine print_formatted_id_bssn_variables(this, points, namefile)

    Interface to print_formatted_id_tpo_variables

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout) :: this

    bssn object to which this PROCEDURE is bound

    double precision, intent(in), optional, DIMENSION(:,:,:,:) :: points
    character(len=*), intent(inout), optional :: namefile

procedure, public, NON_OVERRIDABLE :: print_summary

Prints a summary about the features of the refined mesh

  • interface

    public module subroutine print_summary(this)

    Prints a summary of the properties of the refined mesh, and optionally, to a formatted file whose name is given as the optional argument filename

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this

    Name of the formatted file to print the summary to

procedure, public :: read_bssn_dump_print_formatted

Reads the binary \(\mathrm{ID}\) file printed by compute_and_print_tpo_variables

  • interface

    public module subroutine read_bssn_dump_print_formatted(this, namefile_bin, namefile)

    Interface to read_bssn_dump_print_formatted

    Arguments

    Type IntentOptional Attributes Name
    class(bssn), intent(inout) :: this

    bssn object to which this PROCEDURE is bound

    character(len=*), intent(inout), optional :: namefile_bin
    character(len=*), intent(inout), optional :: namefile

procedure, public, NON_OVERRIDABLE :: setup_standard_tpo_variables

Set up the refined mesh by reading the gravity_grid_parameter.dat parameter file, and read the standard 3+1 \(\mathrm{ID}\) using the given idbase object, on the refined mesh

  • interface

    public module subroutine setup_standard_tpo_variables(tpof, id, dx, dy, dz)

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: tpof
    class(idbase), intent(inout) :: id
    double precision, intent(in), optional :: dx
    double precision, intent(in), optional :: dy
    double precision, intent(in), optional :: dz

procedure, public, NON_OVERRIDABLE :: test_recovery_m2p

Computes the conserved variables from the physical ones, and vice versa, to test that the recovered physical variables are the same to those computed from the \(\mathrm{ID}\). Uses the metric mapped from the grid to the particles.

  • interface

    public module subroutine test_recovery_m2p(this, parts, namefile)

    Tests the recovery. Computes the conserved variables from the physical ones, and then the physical ones from the conserved ones. It then compares the variables computed with the recovery PROCEDURES, with those computed with \(\texttt{SPHINCS_ID}\). Uses the mesh-2-particle.

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(in) :: this

    tpo object which this PROCEDURE is a member of

    class(particles), intent(in) :: parts

    particles object used to map the metric from the mesh to the particles, and to call the recovery procedures

    character(len=*), intent(inout) :: namefile

    Name of the formatted file where the data is printed

Source Code

  TYPE, EXTENDS(tpo):: bssn
  !# TYPE representing the |bssn| formulation of the |ee|


    INTEGER:: call_flag= 0
    !# Flag set to a value different than 0 if the SUBROUTINE
    !  compute_and_print_bssn_variables is called

    !
    !-- Arrays storing the BSSN variables for the |id| on the grid
    !

    TYPE(grid_function):: Gamma_u
    !! Conformal connection \(\bar{\Gamma} ^i_{jk}\)

    TYPE(grid_function_scalar):: phi
    !! Conformal factor \(\phi \)

    TYPE(grid_function_scalar):: trK
    !! Trace of extrinsic curvature \(K \)

    TYPE(grid_function):: A_BSSN3_ll
    !! Conformal traceless extrinsic curvature \(A_{ij} \)

    TYPE(grid_function):: g_BSSN3_ll
    !! Conformal spatial metric \(\gamma_{ij} \)

    TYPE(grid_function):: Ricci_ll
    !! Ricci tensor \(R_{ij} \)

    TYPE(grid_function_scalar):: Ricci_scalar
    !! Ricci scalar \(R^\mu{}_nu\)

    !
    !-- Connection constraints and its l2 norm and loo norm
    !

    TYPE(grid_function):: GC
    !! Connection constraint computed with the |id| on the mesh
    TYPE(grid_function):: GC_parts
    !# Connection constraint computed with the |bssn| |id| on the mesh, and
    !  the hydrodynamical |id| mapped from the particles to the mesh

    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: GC_l2
    !# \(\ell_2\) norm of the connection constraint computed
    !  with the |id| on the mesh
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: GC_parts_l2
    !# \(\ell_2\) norm of the connection constraint computed with the |bssn|
    !  |id| on the mesh, and the hydrodynamical |id| mapped from the particles
    !  to the mesh
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: GC_loo
    !# \(\ell_\infty\) norm, i.e., supremum of the absolute value, of the
    !  connection constraint computed with the |id| on the mesh
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: GC_parts_loo
    !# \(\ell_\infty\) norm, i.e., supremum of the absolute value, of the
    !  connection constraint computed with the |bssn| |id| on the mesh,
    !  and the hydrodynamical |id| mapped from the particles to the mesh
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: GC_int
    !# Integral of the connection constraint computed with the |id| on the mesh
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: GC_parts_int
    !# Integral of the connection constraint computed with the |bssn| |id| on
    !  the mesh, and the hydrodynamical |id| mapped from the particles to the
    !  mesh with the |id| on the mesh

    LOGICAL, PUBLIC:: export_bin
    !# `.TRUE.` if the binary files for SPHINCS_BSSN are to be exported,
    !  `.FALSE.` otherwise
    LOGICAL, PUBLIC:: export_form_xy
    !# `.TRUE.` if the |id| in the formatted files is to be on the xy plane
    !  only, `.FALSE.` otherwise
    LOGICAL, PUBLIC:: export_form_x
    !# `.TRUE.` if the |id| in the formatted files is to be on the x axis
    !  only, `.FALSE.` otherwise

    TYPE(timer):: bssn_computer_timer
    !# Timer that times how long it takes to compute the |bssn| variables on
    !  the refined mesh


    CONTAINS


    !-------------------!
    !--  SUBROUTINES  --!
    !-------------------!


    PROCEDURE:: define_allocate_fields => allocate_bssn_fields
    !! Allocates memory for the [[bssn]] member grid functions

    PROCEDURE:: deallocate_fields => deallocate_bssn_fields
    !! Deallocates memory for the [[bssn]] member arrays

    PROCEDURE:: compute_and_print_tpo_variables &
                    => compute_and_print_bssn_variables
    !# Computes the |bssn| variables at the particle positions, and optionally
    !  prints them to a binary file to be read by \(\texttt{SPHINCS_BSSN}\)
    !  and \(\texttt{splash}\), and to a formatted file to be read by
    !  \(\texttt{gnuplot}\), by calling
    !  [[bssn:print_formatted_id_tpo_variables]]

    PROCEDURE:: print_formatted_id_tpo_variables &
                    => print_formatted_id_bssn_variables
    !! Prints the |bssn| |id| to a formatted file

    PROCEDURE:: compute_and_print_tpo_constraints_grid &
                    => compute_and_print_bssn_constraints_grid
    !# Computes the |bssn| constraints using the full |id| on the refined mesh,
    !  prints a summary with the statistics for the constraints. Optionally,
    !  prints the constraints to a formatted file to be read by
    !  \(\texttt{gnuplot}\), and analyze the constraints by calling
    !  [[tpo:analyze_constraint]]

    PROCEDURE:: compute_and_print_tpo_constraints_particles &
                    => compute_and_print_bssn_constraints_particles
    !# Computes the |bssn| constraints using the |bssn| |id| on the refined
    !  mesh and the hydrodynamical |id| mapped from the particles to the mesh,
    !  prints a summary with the statistics for the constraints. Optionally,
    !  prints the constraints to a formatted file to be read by
    !  \(\texttt{gnuplot}\), and analyze the constraints by calling
    !  [[tpo:analyze_constraint]]


    PROCEDURE:: destruct_bssn
    !# Destructor for the EXTENDED TYPE bssn, not ABSTRACT TYPE tpo

    FINAL    :: destructor
    !# Destructor; finalizes members from both TYPES tpo and bssn,
    !  by calling [[bssn:destruct_bssn]]

    PROCEDURE, PUBLIC:: read_bssn_dump_print_formatted
    !# Reads the binary |id| file printed by
    !  [[bssn:compute_and_print_tpo_variables]]

    PROCEDURE, PUBLIC:: compute_ricci
    !# Computes the Ricci tensor and the Ricci scalar on the mesh


  END TYPE bssn