tpo Derived Type

type, public, ABSTRACT :: tpo

ABSTRACT TYPE representing the standard formulation of the Einstein equations


Inherits

type~~tpo~~InheritsGraph type~tpo tpo grid_function grid_function type~tpo->grid_function coords, shift_u, g_phys3_ll, K_phys3_ll, S, MC, MC_parts, S_parts grid_function_scalar grid_function_scalar type~tpo->grid_function_scalar rad_coord, lapse, HC, HC_parts, rho, rho_parts level level type~tpo->level levels timer timer type~tpo->timer grid_timer, importer_timer

Inherited by

type~~tpo~~InheritedByGraph type~tpo tpo type~bssn bssn type~bssn->type~tpo

Contents

Source Code

tpo

Components

Type Visibility Attributes Name Initial
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 :: 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

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_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

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 :: 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-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(compute_and_print_tpo_constraints_grid_interface), public, deferred :: compute_and_print_tpo_constraints_grid

Computes the constraints specific to the formulation identified by an EXTENDED TYPE, using the full \(\mathrm{ID}\) on the refined mesh

  • subroutine compute_and_print_tpo_constraints_grid_interface(this, id, namefile, name_logfile, points) Prototype

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout), TARGET :: this
    class(idbase), intent(inout) :: id
    character(len=*), intent(inout) :: namefile
    character(len=*), intent(inout) :: name_logfile
    double precision, intent(in), optional, DIMENSION(:,:,:,:), TARGET :: points

procedure(compute_and_print_tpo_constraints_particles_interface), public, deferred :: compute_and_print_tpo_constraints_particles

Computes the constraints specific to the formulation identified by an EXTENDED TYPE, using the \(\mathrm{BSSNOK}\) \(\mathrm{ID}\) on the refined mesh and the hydrodynamical \(\mathrm{ID}\) mapped from the particles to the mesh

  • subroutine compute_and_print_tpo_constraints_particles_interface(this, parts_obj, namefile, name_logfile, points) Prototype

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout), TARGET :: this
    class(particles), intent(inout) :: parts_obj
    character(len=*), intent(inout) :: namefile
    character(len=*), intent(inout) :: name_logfile
    double precision, intent(in), optional, DIMENSION(:,:,:,:), TARGET :: points

procedure(compute_and_print_tpo_variables_interface), public, deferred :: compute_and_print_tpo_variables

Compute the fields specific to the formulation identified by an EXTENDED TYPE, starting from the standard 3+1 fields

  • subroutine compute_and_print_tpo_variables_interface(this, namefile) Prototype

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    character(len=*), intent(inout), optional :: namefile

procedure(deallocate_fields_interface), public, deferred :: deallocate_fields

Deallocates memory for the fields specific to the formulation identified by an EXTENDED TYPE

  • subroutine deallocate_fields_interface(this) Prototype

    Arguments

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

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(define_allocate_fields_interface), public, deferred :: define_allocate_fields

Allocates memory for the fields specific to the formulation identified by an EXTENDED TYPE

  • subroutine define_allocate_fields_interface(this) Prototype

    Arguments

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

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(print_formatted_id_tpo_variables_interface), public, deferred :: print_formatted_id_tpo_variables

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

  • subroutine print_formatted_id_tpo_variables_interface(this, points, namefile) Prototype

    Arguments

    Type IntentOptional Attributes Name
    class(tpo), intent(inout) :: this
    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, 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, ABSTRACT:: tpo
  !! ABSTRACT TYPE representing the standard \(3+1\) formulation of the |ee|

    INTEGER:: tpo_id_number
    !! Negative integer that identifies the [[tpo]] object

    INTEGER:: n_matter
    !! Number of matter objects in the physical system

    INTEGER, 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

    !
    !-- Mesh variables
    !

    INTEGER:: nlevels
    !! Number of refinement levels

    TYPE(level), DIMENSION(:), ALLOCATABLE :: levels
    !! Array containing the information on each refinement level
    TYPE(grid_function):: coords
    !! Grid function storing the Cartesian coordinates
    TYPE(grid_function_scalar):: rad_coord
    !! Grid scalar function storing the radial coordinates of each grid point

    !
    !-- ADM fields
    !

    TYPE(grid_function_scalar):: lapse
    !! Grid scalar function storing the lapse function
    TYPE(grid_function):: shift_u
    !! Grid function storing the shift vector
    TYPE(grid_function):: g_phys3_ll
    !! Grid function storing the spatial metric
    TYPE(grid_function):: K_phys3_ll
    !! Grid function storing the extrinsic curvature

    !
    !-- Constraint violations
    !

    TYPE(grid_function_scalar):: HC
    !# Grid scalar function storing the Hamiltonian constraint (violations)
    !  computed using the ID on the mesh
    TYPE(grid_function_scalar):: HC_parts
    !# Grid scalar function storing the Hamiltonian constraint (violations)
    !  computed using the stress-energy tensor mapped from the particles to the
    !  mesh
    TYPE(grid_function_scalar):: rho
    !# Grid scalar function storing the matter source in the Hamiltonian
    !  constraint computed using the ID on the mesh, multiplied by \(16\pi\)
    TYPE(grid_function):: S
    !# Grid function storing the matter source in the momentum
    !  constraint computed using the ID on the mesh, multiplied by \(8\pi\)
    TYPE(grid_function):: MC
    !# Grid function storing the momentum constraint (violations)
    !  computed using the ID on the mesh
    TYPE(grid_function):: MC_parts
    !# Grid function storing the momentum constraint (violations)
    !  computed using the stress-energy tensor mapped from the particles to the
    !  mesh
    TYPE(grid_function_scalar):: 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 \(16\pi\)
    TYPE(grid_function):: 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 \(8\pi\)

    !
    !-- Norms of constraint violations
    !

    DOUBLE PRECISION, DIMENSION(:),   ALLOCATABLE:: HC_l2
    !! \(\ell_2\) norm of the Hamiltonian constraint computed on the mesh
    DOUBLE PRECISION, DIMENSION(:),   ALLOCATABLE:: HC_parts_l2
    !# \(\ell_2\) norm of the Hamiltonian constraint computed on the mesh, using
    !  the stress-energy tensor mapped from the particles
    DOUBLE PRECISION, DIMENSION(:),   ALLOCATABLE:: HC_loo
    !# \(\ell_\infty\) norm of the Hamiltonian constraint computed on the mesh
    !  (i.e., its maximum)
    DOUBLE PRECISION, DIMENSION(:),   ALLOCATABLE:: HC_parts_loo
    !# \(\ell_\infty\) norm of the Hamiltonian constraint computed on the mesh,
    !  using the stress-energy tensor mapped from the particles
    !  (i.e., its maximum)
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: MC_l2
    !! \(\ell_2\) norm of the momentum constraint computed on the mesh
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: MC_parts_l2
    !# \(\ell_2\) norm of the Hamiltonian constraint computed on the mesh, using
    !  the stress-energy tensor mapped from the particles
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: MC_loo
    !# \(\ell_\infty\) norm of the momentum constraint computed on the mesh
    !  (i.e., its maximum)
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: MC_parts_loo
    !# \(\ell_\infty\) norm of the momentum constraint computed on the mesh,
    !  using the stress-energy tensor mapped from the particles
    !  (i.e., its maximum)
    DOUBLE PRECISION, DIMENSION(:),   ALLOCATABLE:: HC_int
    !# Integral the Hamiltonian constraint computed using the stress-energy
    !  tensor mapped from the particles
    DOUBLE PRECISION, DIMENSION(:),   ALLOCATABLE:: HC_parts_int
    !! Integral the Hamiltonian constraint computed on the mesh
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: MC_int
    !! Integral of the momentum constraint computed on the mesh
    DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE:: MC_parts_int
    !# Integral of the momentum constraint computed using the stress-energy
    !  tensor mapped from the particles

    !
    !-- Steering variables
    !

    ! Variables to decide if and how to export the constraints

    INTEGER, PUBLIC:: cons_step
    !! Constraint violations are printed to file every cons_step-th grid point
    !  along each Cartesian direction
    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 \([0,10^{-7}],[10^3,\infty],[10^n,10^m]\), with
    !  \(n\in\{-7,2\}\) and \(m=n+1\), are to be printed to file;
    !  `.FALSE.` otherwise
    LOGICAL, PUBLIC:: export_constraints_xy
    !! `.TRUE.` if only the constrain violations on the \(xy\) plane are to be
    !  printed to file, `.FALSE.` otherwise
    LOGICAL, PUBLIC:: export_constraints_x
    !! `.TRUE.` if only the constrain violations on the \(x\) axis are to be
    !  printed to file, `.FALSE.` otherwise

    !
    !-- Timers
    !

    TYPE(timer):: grid_timer
    !# Timer that times how long it takes to set up the grid and allocate
    !  the grid functions
    TYPE(timer):: importer_timer
    !# Timer that times how long it takes to import the |lorene| ID
    !  on the mesh


    CONTAINS


    !-----------------------------------!
    !--  NON_OVERRIDABLE SUBROUTINES  --!
    !-----------------------------------!


    PROCEDURE, 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 |id| using the given
    !  [[idbase]] object, on the refined mesh

    PROCEDURE, NON_OVERRIDABLE:: deallocate_standard_tpo_variables
    !! Deallocates memory for the standard 3+1 fields

    PROCEDURE, 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
    !  @todo complete this documentation entries with more details

    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 |id|. Uses the metric mapped from the grid to the
    !  particles. @todo add reference for recovery

    PROCEDURE, PUBLIC, NON_OVERRIDABLE:: compute_adm_momentum_fluid_m2p
    !# Computes an estimate of the \(\mathrm{ADM}\) linear momentum of the
    !  fluid using the |sph| hydro fields, and the spacetime metric mapped
    !  from the mesh
    !  @todo add reference

    PROCEDURE, NON_OVERRIDABLE:: print_summary
    !# Prints a summary about the features of the refined mesh


    !----------------------------!
    !--  DEFERRED SUBROUTINES  --!
    !----------------------------!


    PROCEDURE(define_allocate_fields_interface), &
      DEFERRED:: define_allocate_fields
    !# Allocates memory for the fields specific to the formulation identified
    !  by an EXTENDED TYPE

    PROCEDURE(deallocate_fields_interface), DEFERRED:: deallocate_fields
    !# Deallocates memory for the fields specific to the formulation identified
    !  by an EXTENDED TYPE

    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 |id|
    !  on the refined mesh, or the spacetime |id| on the refined mesh and
    !  the hydrodynamical |id| mapped from the particles to the refined mesh

    PROCEDURE(compute_and_print_tpo_constraints_grid_interface), &
      DEFERRED:: compute_and_print_tpo_constraints_grid
    !# Computes the constraints specific to the formulation identified by an
    !  EXTENDED TYPE, using the full |id| on the refined mesh

    PROCEDURE(compute_and_print_tpo_constraints_particles_interface), &
      DEFERRED:: compute_and_print_tpo_constraints_particles
    !# Computes the constraints specific to the formulation identified by an
    !  EXTENDED TYPE, using the |bssn| |id| on the refined
    !  mesh and the hydrodynamical |id| mapped from the particles to the mesh

    PROCEDURE(compute_and_print_tpo_variables_interface), PUBLIC, &
      DEFERRED:: compute_and_print_tpo_variables
    !# Compute the fields specific to the formulation identified by an
    !  EXTENDED TYPE, starting from the standard 3+1 fields

    PROCEDURE(print_formatted_id_tpo_variables_interface), PUBLIC, &
      DEFERRED:: print_formatted_id_tpo_variables
    !! Prints the spacetime |id| to a formatted file


    !-----------------!
    !--  FUNCTIONS  --!
    !-----------------!


    PROCEDURE:: abs_values_in

    PROCEDURE, PUBLIC:: get_grid_point

    PROCEDURE, PUBLIC:: get_nlevels

    PROCEDURE, PUBLIC:: get_levels

    PROCEDURE, PUBLIC:: get_dx

    PROCEDURE, PUBLIC:: get_dy

    PROCEDURE, PUBLIC:: get_dz

    PROCEDURE, PUBLIC:: get_ngrid_x

    PROCEDURE, PUBLIC:: get_ngrid_y

    PROCEDURE, PUBLIC:: get_ngrid_z

    PROCEDURE, PUBLIC:: get_xR

    PROCEDURE, PUBLIC:: get_yR

    PROCEDURE, PUBLIC:: get_zR

    PROCEDURE, PUBLIC:: get_HC

    PROCEDURE, PUBLIC:: get_MC

    PROCEDURE, PUBLIC:: get_HC_parts

    PROCEDURE, PUBLIC:: get_MC_parts


  END TYPE tpo