Bimetric 3+1 toolkit for spherical symmetry
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
Implementation Details

Main Class

  • BimetricEvolve: Encapsulates a 3+1 evolution solver for bimetric spacetimes.


Modules across the source files

     +-------------------------------------------------------+    +-----------------+
     |                  Bimetric Evolution                   |    |  Bimetric Model |
     +-------------------------------------------------------+    +-----------------+
     | bim-solver.cpp   eomEvolution.h    eomGauge.h         |    | bimetricModel.h |
     | maximalSlice.h   eomMatter.h       eomLapseRatios.h   |    +-----------------+
     |                  eomSources.h      eomMiscEquations.h |
     |                  eomConstraints.h                     |
     +-------------------------------------------------------+
   
  +-----------------+   +-------------------+   +----------------+   +---------------+
  |   Grid Driver   |   |   Initial Data    |   | MoL Integrator |   |  Output Data  |
  +-----------------+   +-------------------+   +----------------+   +---------------+
  | gridDriver.h    |   | gridInitialData.h |   |  integrator.h  |   | gridOutput.h  |
  | gridFunctions.h |   +-------------------+   +----------------+   +---------------+
  +-----------------+
         +---------------------+      +--------------+      +-----------------+
         |  Numerical Methods  |      |  Data Types  |      |    Utilities    |
         +---------------------+      +--------------+      +-----------------+
         | finiteDifferences.h |      | matrix.h     |      | trackUsedTime.h |
         | methodOfLines.h     |      | dataTypes.h  |      | mpiWorld.h      |
         | cubicSpline.h       |      +--------------+      | mpiDummyWorld.h |
         | bandSol.h           |                            | paramsHolder.h  |
         +---------------------+                            | slog.h          |
                                                            +-----------------+ 

Data Flow

  +---------------------------------+              +----------------------------+
  | Mathematica Notebook            |------------>>| Initial Data & Parameters  |
  | bimetric-ss-in-3+1 (cpp solver) |  Solve the   |  input.dat  (initial data) |
  +---------------------------------+  constraint  |  config.ini (parameters)   |
            |   Export the             equations   +----------------------------+
            |   C++ code                                         |
            V                                                    V
  +---------------------------------+  eom-std.h   +----------------------------+
  | eomEvolution.h eomMatter.h      |------------>>|       bim-solver.cpp       |
  | eomSources.h   eomConstraints.h |              +----------------------------+
  | eomGauge.h     eomLapseRatios.h |                            |
  | eomMiscEquations.h              |                            V
  +---------------------------------+              +----------------------------+
                                                   |         output.dat         |
     See also:  maximalSlice.h                     +----------------------------+
      
     Numerics:  methodOfLines.h, finiteDifferences.h,
                dataTypes.h, matrix.h, bandSol.h, cubicSpline.h  

Regularization of spherically symmetric evolution codes

  • A nice overview of the parity conditions at r = 0 can be found in [11].
  • Also, parity conditions for components of vectors and tensors can be found in Table 1 in [4].


Finite Difference Approximation

At the core of finite difference approximation is a discretization of the spacetime, or a numerical grid (shown below). Then, a function f(t,r) is represented by values at a discrete set of points (m,n).

                p in M          <-->     (t,r)      <-->        (m,n)
         (point in a manifold)        (in a chart)        (discretized point)
                                           |                      |
                                           V                      V
                                         f(t,r)                 f(m,n)    
The points (m,n) are called the grid points or nodes. The grid points are chosen to reside at the center of grid cells (the other approach would be to associate the grid points with the vertices of the grid cells). The data structure GridPoint is a wrapper for the set of fields {f,...} that arekept in the memory for the point (m,n).The grid consists of the interior points, and the virtual points at the boundary. The virtual points are used to impose the boundary conditions. A cell-centered uniform grid with nLen grid points in each grid row is shown below (see also [3], p.192).
         left virtual                grid cell        grid point     right virtual
            point(s)                 .---^---.            |             point(s)
      .        |               n-1  '    n    '  n+1      |                |
      :   +----|----+-------+-------+---------+-------+---|---+-------+----|----+
     m+1  |    O    |       |       | (m+1,n) |       |   X   |       |    O    |
          +---------+-------+-------.=========.-------+-------+-------+---------+
      m   |         |       |(m,n-1)|  (m,n)  |(m,n+1)|       |       |         |
          +---------+-------+-------'========='-------+-------+-------+---------+
     m-1  |         |       |       |         |       |       |       |         |
      .   +---------'-------+-------+---------+-------+-------+-------'---------+
      :        0    |   1       2      . . .                    nLen  |  nLen+1
                    ^                                                 ^
                   r_min                                            r_max       
The cell-centered uniform grid implies:

     r = r_min + ( n - 1/2 ) delta_r
     delta_r = ( r_max - r_min ) / nLen
  • if nGhost = 1, n ranges between 1 and nLen; the boundary virtual points are grid[m][0] and grid[m][nLen+1]
  • m ranges between 1 and mLen, and wraps around; the slices grid[mLen+1]... are used for intermediate integration steps.

In parallel simulations, the grid is split across processes using MPI, allowing storage of large grid functions and parallel work on each part.