Notebook Status: Validated
Validation Notes: This module has been validated to exhibit convergence to zero of the Hamiltonian constraint violation at the expected order to the exact solution after a short numerical evolution of the initial data (see plots at bottom), and all quantities have been validated against the original SENR code. Finally, excellent agreement is seen in the gravitational wave signal from the ringing remnant black hole for multiple decades in amplitude when compared to black hole perturbation theory predictions.
Here we use NRPy+ to generate the C source code necessary to set up initial data for two black holes (Brill-Lindquist, Brill & Lindquist, Phys. Rev. 131, 471, 1963; see also Eq. 1 of Brandt & Brügmann, arXiv:gr-qc/9711015v1). Then we use it to generate the RHS expressions for Method of Lines time integration based on the explicit Runge-Kutta fourth-order scheme (RK4).
The entire algorithm is outlined as follows, with links to the relevant NRPy+ tutorial notebooks listed at each step:
This notebook is organized as follows
BSSN.BrillLindquist
NRPy+ modulel_max
(set to 2 here)free_parameters.h
BrillLindquist_Playground.c
: The Main C Code# Step P0: Set the option for evolving the initial data forward in time
# Options include "low resolution" and "high resolution"
EvolOption = "low resolution"
# Step P1: Import needed NRPy+ core modules:
from outputC import lhrh,outCfunction,outC_function_dict # NRPy+: Core C code output module
import finite_difference as fin # NRPy+: Finite difference C code generation module
import NRPy_param_funcs as par # NRPy+: Parameter interface
import grid as gri # NRPy+: Functions having to do with numerical grids
import indexedexp as ixp # NRPy+: Symbolic indexed expression (e.g., tensors, vectors, etc.) support
import reference_metric as rfm # NRPy+: Reference metric support
import deprecated_reference_metric as evil_rfm # NRPy+: PLEASE DON'T USE.
import cmdline_helper as cmd # NRPy+: Multi-platform Python command-line interface
import shutil, os, sys, time # Standard Python modules for multiplatform OS-level functions, benchmarking
# Step P2: Create C code output directory:
Ccodesdir = os.path.join("BSSN_Two_BHs_Collide_Ccodes_psi4")
# First remove C code output directory if it exists
# Courtesy https://stackoverflow.com/questions/303200/how-do-i-remove-delete-a-folder-that-is-not-empty
shutil.rmtree(Ccodesdir, ignore_errors=True)
# Then create a fresh directory
cmd.mkdir(Ccodesdir)
# Step P3: Create executable output directory:
outdir = os.path.join(Ccodesdir, "output")
cmd.mkdir(outdir)
# Step 1: Set the spatial dimension parameter
# to three (BSSN is a 3+1 decomposition
# of Einstein's equations), and then read
# the parameter as DIM.
DIM = 3
par.set_parval_from_str("grid::DIM",DIM)
# Step 1.a: Enable SIMD-optimized code?
# I.e., generate BSSN and Ricci C code kernels using SIMD-vectorized
# compiler intrinsics, which *greatly improve the code's performance*,
# though at the expense of making the C-code kernels less
# human-readable.
# * Important note in case you wish to modify the BSSN/Ricci kernels
# here by adding expressions containing transcendental functions
# (e.g., certain scalar fields):
# Note that SIMD-based transcendental function intrinsics are not
# supported by the default installation of gcc or clang (you will
# need to use e.g., the SLEEF library from sleef.org, for this
# purpose). The Intel compiler suite does support these intrinsics
# however without the need for external libraries.
enable_SIMD = True
# Step 2: Set some core parameters, including CoordSystem MoL timestepping algorithm,
# FD order, floating point precision, and CFL factor:
# Choices are: Spherical, SinhSpherical, SinhSphericalv2, Cylindrical, SinhCylindrical,
# SymTP, SinhSymTP
CoordSystem = "SinhSpherical"
# Decompose psi_4 (second time derivative of gravitational
# wave strain) into all spin-weight=-2
# l,m spherical harmonics, starting at l=2
# going up to and including l_max, set here:
l_max = 2
if EvolOption == "low resolution":
# domain_size sets the default value for:
# * Spherical's params.RMAX
# * SinhSpherical*'s params.AMAX
# * Cartesians*'s -params.{x,y,z}min & .{x,y,z}max
# * Cylindrical's -params.ZMIN & .{Z,RHO}MAX
# * SinhCylindrical's params.AMPL{RHO,Z}
# * *SymTP's params.AMAX
domain_size = 150 # Length scale of computational domain
FD_order = 8 # Finite difference order: even numbers only, starting with 2. 12 is generally unstable
elif EvolOption == "high resolution":
# See above for description of the domain_size parameter
domain_size = 300 # Length scale of computational domain
FD_order = 10 # Finite difference order: even numbers only, starting with 2. 12 is generally unstable
else:
print("Error: EvolOption = "+str(EvolOption)+" unrecognized.")
exit(1)
# sinh_width sets the default value for:
# * SinhSpherical's params.SINHW
# * SinhCylindrical's params.SINHW{RHO,Z}
# * SinhSymTP's params.SINHWAA
sinh_width = 0.2 # If Sinh* coordinates chosen
# sinhv2_const_dr sets the default value for:
# * SinhSphericalv2's params.const_dr
# * SinhCylindricalv2's params.const_d{rho,z}
sinhv2_const_dr = 0.05 # If Sinh*v2 coordinates chosen
# SymTP_bScale sets the default value for:
# * SinhSymTP's params.bScale
SymTP_bScale = 0.5 # If SymTP chosen
# Step 2.b: Set the timestepping order,
# the core data type, and the CFL factor.
# Step 2.b: Set the order of spatial and temporal derivatives;
# the core data type, and the CFL factor.
# RK_method choices include: Euler, "RK2 Heun", "RK2 MP", "RK2 Ralston", RK3, "RK3 Heun", "RK3 Ralston",
# SSPRK3, RK4, DP5, DP5alt, CK5, DP6, L6, DP8
RK_method = "RK4"
REAL = "double" # Best to use double here.
CFL_FACTOR= 1.0 # (GETS OVERWRITTEN WHEN EXECUTED.) In pure axisymmetry (symmetry_axes = 2 below) 1.0 works fine. Otherwise 0.5 or lower.
# Step 3: Generate Runge-Kutta-based (RK-based) timestepping code.
# As described above the Table of Contents, this is a 3-step process:
# 3.A: Evaluate RHSs (RHS_string)
# 3.B: Apply boundary conditions (post_RHS_string, pt 1)
# 3.C: Enforce det(gammabar) = det(gammahat) constraint (post_RHS_string, pt 2)
import MoLtimestepping.C_Code_Generation as MoL
from MoLtimestepping.RK_Butcher_Table_Dictionary import Butcher_dict
RK_order = Butcher_dict[RK_method][1]
cmd.mkdir(os.path.join(Ccodesdir,"MoLtimestepping/"))
MoL.MoL_C_Code_Generation(RK_method,
RHS_string = """
Ricci_eval(&rfmstruct, ¶ms, RK_INPUT_GFS, auxevol_gfs);
rhs_eval(&rfmstruct, ¶ms, auxevol_gfs, RK_INPUT_GFS, RK_OUTPUT_GFS);""",
post_RHS_string = """
apply_bcs_curvilinear(¶ms, &bcstruct, NUM_EVOL_GFS, evol_gf_parity, RK_OUTPUT_GFS);
enforce_detgammahat_constraint(&rfmstruct, ¶ms, RK_OUTPUT_GFS);\n""",
outdir = os.path.join(Ccodesdir,"MoLtimestepping/"))
# Step 4: Set the coordinate system for the numerical grid
par.set_parval_from_str("reference_metric::CoordSystem",CoordSystem)
rfm.reference_metric() # Create ReU, ReDD needed for rescaling B-L initial data, generating BSSN RHSs, etc.
# Step 5: Set the finite differencing order to FD_order (set above).
par.set_parval_from_str("finite_difference::FD_CENTDERIVS_ORDER", FD_order)
# Step 6: Copy SIMD/SIMD_intrinsics.h to $Ccodesdir/SIMD/SIMD_intrinsics.h
cmd.mkdir(os.path.join(Ccodesdir,"SIMD"))
shutil.copy(os.path.join("SIMD/")+"SIMD_intrinsics.h",os.path.join(Ccodesdir,"SIMD/"))
# Step 7: Set the direction=2 (phi) axis to be the symmetry axis; i.e.,
# axis "2", corresponding to the i2 direction.
# This sets all spatial derivatives in the phi direction to zero.
par.set_parval_from_str("indexedexp::symmetry_axes","2")
In order for our explicit-timestepping numerical solution to the scalar wave equation to be stable, it must satisfy the CFL condition: $$ \Delta t \le \frac{\min(ds_i)}{c}, $$ where $c$ is the wavespeed, and $$ds_i = h_i \Delta x^i$$ is the proper distance between neighboring gridpoints in the $i$th direction (in 3D, there are 3 directions), $h_i$ is the $i$th reference metric scale factor, and $\Delta x^i$ is the uniform grid spacing in the $i$th direction:
# Output the find_timestep() function to a C file.
import deprecated_reference_metric as evil_rfm
evil_rfm.out_timestep_func_to_file(os.path.join(Ccodesdir,"find_timestep.h"))
BSSN.BrillLindquist
NRPy+ module [Back to top]¶The BSSN.BrillLindquist
NRPy+ module does the following:
import BSSN.BrillLindquist as bl
def BrillLindquistID():
print("Generating optimized C code for Brill-Lindquist initial data. May take a while, depending on CoordSystem.")
start = time.time()
bl.BrillLindquist() # Registers ID C function in dictionary, used below to output to file.
with open(os.path.join(Ccodesdir,"initial_data.h"),"w") as file:
file.write(outC_function_dict["initial_data"])
end = time.time()
print("(BENCH) Finished BL initial data codegen in "+str(end-start)+" seconds.")
import BSSN.BSSN_RHSs as rhs
import BSSN.BSSN_gauge_RHSs as gaugerhs
# Set the *covariant*, second-order Gamma-driving shift condition
par.set_parval_from_str("BSSN.BSSN_gauge_RHSs::ShiftEvolutionOption", "GammaDriving2ndOrder_Covariant")
print("Generating symbolic expressions for BSSN RHSs...")
start = time.time()
# Enable rfm_precompute infrastructure, which results in
# BSSN RHSs that are free of transcendental functions,
# even in curvilinear coordinates, so long as
# ConformalFactor is set to "W" (default).
cmd.mkdir(os.path.join(Ccodesdir,"rfm_files/"))
par.set_parval_from_str("reference_metric::enable_rfm_precompute","True")
par.set_parval_from_str("reference_metric::rfm_precompute_Ccode_outdir",os.path.join(Ccodesdir,"rfm_files/"))
# Evaluate BSSN + BSSN gauge RHSs with rfm_precompute enabled:
import BSSN.BSSN_quantities as Bq
par.set_parval_from_str("BSSN.BSSN_quantities::LeaveRicciSymbolic","True")
rhs.BSSN_RHSs()
gaugerhs.BSSN_gauge_RHSs()
# We use betaU as our upwinding control vector:
Bq.BSSN_basic_tensors()
betaU = Bq.betaU
import BSSN.Enforce_Detgammahat_Constraint as EGC
enforce_detg_constraint_symb_expressions = EGC.Enforce_Detgammahat_Constraint_symb_expressions()
# Next compute Ricci tensor
par.set_parval_from_str("BSSN.BSSN_quantities::LeaveRicciSymbolic","False")
Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()
# Now register the Hamiltonian as a gridfunction.
H = gri.register_gridfunctions("AUX","H")
# Then define the Hamiltonian constraint and output the optimized C code.
import BSSN.BSSN_constraints as bssncon
bssncon.BSSN_constraints(add_T4UUmunu_source_terms=False)
# Now that we are finished with all the rfm hatted
# quantities in generic precomputed functional
# form, let's restore them to their closed-
# form expressions.
par.set_parval_from_str("reference_metric::enable_rfm_precompute","False") # Reset to False to disable rfm_precompute.
rfm.ref_metric__hatted_quantities()
end = time.time()
print("(BENCH) Finished BSSN symbolic expressions in "+str(end-start)+" seconds.")
def BSSN_RHSs():
print("Generating C code for BSSN RHSs in "+par.parval_from_str("reference_metric::CoordSystem")+" coordinates.")
start = time.time()
# Construct the left-hand sides and right-hand-side expressions for all BSSN RHSs
lhs_names = [ "alpha", "cf", "trK"]
rhs_exprs = [gaugerhs.alpha_rhs, rhs.cf_rhs, rhs.trK_rhs]
for i in range(3):
lhs_names.append( "betU"+str(i))
rhs_exprs.append(gaugerhs.bet_rhsU[i])
lhs_names.append( "lambdaU"+str(i))
rhs_exprs.append(rhs.lambda_rhsU[i])
lhs_names.append( "vetU"+str(i))
rhs_exprs.append(gaugerhs.vet_rhsU[i])
for j in range(i,3):
lhs_names.append( "aDD"+str(i)+str(j))
rhs_exprs.append(rhs.a_rhsDD[i][j])
lhs_names.append( "hDD"+str(i)+str(j))
rhs_exprs.append(rhs.h_rhsDD[i][j])
# Sort the lhss list alphabetically, and rhss to match.
# This ensures the RHSs are evaluated in the same order
# they're allocated in memory:
lhs_names,rhs_exprs = [list(x) for x in zip(*sorted(zip(lhs_names,rhs_exprs), key=lambda pair: pair[0]))]
# Declare the list of lhrh's
BSSN_evol_rhss = []
for var in range(len(lhs_names)):
BSSN_evol_rhss.append(lhrh(lhs=gri.gfaccess("rhs_gfs",lhs_names[var]),rhs=rhs_exprs[var]))
# Set up the C function for the BSSN RHSs
desc="Evaluate the BSSN RHSs"
name="rhs_eval"
outCfunction(
outfile = os.path.join(Ccodesdir,name+".h"), desc=desc, name=name,
params = """rfm_struct *restrict rfmstruct,const paramstruct *restrict params,
const REAL *restrict auxevol_gfs,const REAL *restrict in_gfs,REAL *restrict rhs_gfs""",
body = fin.FD_outputC("returnstring",BSSN_evol_rhss, params="outCverbose=False,enable_SIMD=True",
upwindcontrolvec=betaU),
loopopts = "InteriorPoints,enable_SIMD,enable_rfm_precompute")
end = time.time()
print("(BENCH) Finished BSSN_RHS C codegen in " + str(end - start) + " seconds.")
def Ricci():
print("Generating C code for Ricci tensor in "+par.parval_from_str("reference_metric::CoordSystem")+" coordinates.")
start = time.time()
desc="Evaluate the Ricci tensor"
name="Ricci_eval"
outCfunction(
outfile = os.path.join(Ccodesdir,name+".h"), desc=desc, name=name,
params = """rfm_struct *restrict rfmstruct,const paramstruct *restrict params,
const REAL *restrict in_gfs,REAL *restrict auxevol_gfs""",
body = fin.FD_outputC("returnstring",
[lhrh(lhs=gri.gfaccess("auxevol_gfs","RbarDD00"),rhs=Bq.RbarDD[0][0]),
lhrh(lhs=gri.gfaccess("auxevol_gfs","RbarDD01"),rhs=Bq.RbarDD[0][1]),
lhrh(lhs=gri.gfaccess("auxevol_gfs","RbarDD02"),rhs=Bq.RbarDD[0][2]),
lhrh(lhs=gri.gfaccess("auxevol_gfs","RbarDD11"),rhs=Bq.RbarDD[1][1]),
lhrh(lhs=gri.gfaccess("auxevol_gfs","RbarDD12"),rhs=Bq.RbarDD[1][2]),
lhrh(lhs=gri.gfaccess("auxevol_gfs","RbarDD22"),rhs=Bq.RbarDD[2][2])],
params="outCverbose=False,enable_SIMD=True"),
loopopts = "InteriorPoints,enable_SIMD,enable_rfm_precompute")
end = time.time()
print("(BENCH) Finished Ricci C codegen in " + str(end - start) + " seconds.")
Generating symbolic expressions for BSSN RHSs... (BENCH) Finished BSSN symbolic expressions in 4.205245494842529 seconds.
Next output the C code for evaluating the Hamiltonian constraint (Tutorial). In the absence of numerical error, this constraint should evaluate to zero. However it does not due to numerical (typically truncation and roundoff) error. We will therefore measure the Hamiltonian constraint violation to gauge the accuracy of our simulation, and, ultimately determine whether errors are dominated by numerical finite differencing (truncation) error as expected.
def Hamiltonian():
start = time.time()
print("Generating optimized C code for Hamiltonian constraint. May take a while, depending on CoordSystem.")
# Set up the C function for the Hamiltonian RHS
desc="Evaluate the Hamiltonian constraint"
name="Hamiltonian_constraint"
outCfunction(
outfile = os.path.join(Ccodesdir,name+".h"), desc=desc, name=name,
params = """rfm_struct *restrict rfmstruct,const paramstruct *restrict params,
REAL *restrict in_gfs, REAL *restrict aux_gfs""",
body = fin.FD_outputC("returnstring",lhrh(lhs=gri.gfaccess("aux_gfs", "H"), rhs=bssncon.H),
params="outCverbose=False"),
loopopts = "InteriorPoints,enable_rfm_precompute")
end = time.time()
print("(BENCH) Finished Hamiltonian C codegen in " + str(end - start) + " seconds.")
Then enforce conformal 3-metric $\det{\bar{\gamma}_{ij}}=\det{\hat{\gamma}_{ij}}$ constraint (Eq. 53 of Ruchlin, Etienne, and Baumgarte (2018)), as documented in the corresponding NRPy+ tutorial notebook
Applying curvilinear boundary conditions should affect the initial data at the outer boundary, and will in general cause the $\det{\bar{\gamma}_{ij}}=\det{\hat{\gamma}_{ij}}$ constraint to be violated there. Thus after we apply these boundary conditions, we must always call the routine for enforcing the $\det{\bar{\gamma}_{ij}}=\det{\hat{\gamma}_{ij}}$ constraint:
def gammadet():
start = time.time()
print("Generating optimized C code for gamma constraint. May take a while, depending on CoordSystem.")
# Set up the C function for the det(gammahat) = det(gammabar)
EGC.output_Enforce_Detgammahat_Constraint_Ccode(Ccodesdir,exprs=enforce_detg_constraint_symb_expressions)
end = time.time()
print("(BENCH) Finished gamma constraint C codegen in " + str(end - start) + " seconds.")
The Weyl scalar $\psi_4$ encodes gravitational wave information in our numerical relativity calculations. For more details on how it is computed, see this NRPy+ tutorial notebook for information on $\psi_4$ and this one on the Quasi-Kinnersley tetrad (as implemented in Baker, Campanelli, Lousto (2001)).
$\psi_4$ is related to the gravitational wave strain via $$ \psi_4 = \ddot{h}_+ - i \ddot{h}_\times, $$ where $\ddot{h}_+$ is the second time derivative of the $+$ polarization of the gravitational wave strain $h$, and $\ddot{h}_\times$ is the second time derivative of the $\times$ polarization of the gravitational wave strain $h$.
import BSSN.Psi4_tetrads as BP4t
par.set_parval_from_str("BSSN.Psi4_tetrads::TetradChoice","QuasiKinnersley")
#par.set_parval_from_str("BSSN.Psi4_tetrads::UseCorrectUnitNormal","True")
import BSSN.Psi4 as BP4
print("Generating symbolic expressions for psi4...")
start = time.time()
BP4.Psi4()
end = time.time()
print("(BENCH) Finished psi4 symbolic expressions in "+str(end-start)+" seconds.")
psi4r_0pt = gri.register_gridfunctions("AUX","psi4r_0pt")
psi4r_1pt = gri.register_gridfunctions("AUX","psi4r_1pt")
psi4r_2pt = gri.register_gridfunctions("AUX","psi4r_2pt")
psi4i_0pt = gri.register_gridfunctions("AUX","psi4i_0pt")
psi4i_1pt = gri.register_gridfunctions("AUX","psi4i_1pt")
psi4i_2pt = gri.register_gridfunctions("AUX","psi4i_2pt")
desc="""Since it's so expensive to compute, instead of evaluating
psi_4 at all interior points, this functions evaluates it on a
point-by-point basis."""
name="psi4"
outCfunction(
outfile = os.path.join(Ccodesdir,name+".h"), desc=desc, name=name,
params = """const paramstruct *restrict params,
const int i0,const int i1,const int i2,
REAL *restrict xx[3], const REAL *restrict in_gfs, REAL *restrict aux_gfs""",
body = """
const int idx = IDX3S(i0,i1,i2);
const REAL xx0 = xx[0][i0];const REAL xx1 = xx[1][i1];const REAL xx2 = xx[2][i2];
// Real part of psi_4, divided into 3 terms
{
#include "Psi4re_pt0_lowlevel.h"
}
{
#include "Psi4re_pt1_lowlevel.h"
}
{
#include "Psi4re_pt2_lowlevel.h"
}
// Imaginary part of psi_4, divided into 3 terms
{
#include "Psi4im_pt0_lowlevel.h"
}
{
#include "Psi4im_pt1_lowlevel.h"
}
{
#include "Psi4im_pt2_lowlevel.h"
}""")
def Psi4re(part):
print("Generating C code for psi4_re_pt"+str(part)+" in "+par.parval_from_str("reference_metric::CoordSystem")+" coordinates.")
start = time.time()
Psi4re_pt = fin.FD_outputC("returnstring",
[lhrh(lhs=gri.gfaccess("aux_gfs","psi4r_"+str(part)+"pt"),rhs=BP4.psi4_re_pt[part])],
params="outCverbose=False,CSE_sorting=none") # Generating the CSE for psi4 is the slowest
# operation in this notebook, and much of the CSE
# time is spent sorting CSE expressions. Disabling
# this sorting makes the C codegen 3-4x faster,
# but the tradeoff is that every time this is
# run, the CSE patterns will be different
# (though they should result in mathematically
# *identical* expressions). You can expect
# roundoff-level differences as a result.
with open(os.path.join(Ccodesdir,"Psi4re_pt"+str(part)+"_lowlevel.h"), "w") as file:
file.write(Psi4re_pt)
end = time.time()
print("(BENCH) Finished generating psi4_re_pt"+str(part)+" in "+str(end-start)+" seconds.")
def Psi4im(part):
print("Generating C code for psi4_im_pt"+str(part)+" in "+par.parval_from_str("reference_metric::CoordSystem")+" coordinates.")
start = time.time()
Psi4im_pt = fin.FD_outputC("returnstring",
[lhrh(lhs=gri.gfaccess("aux_gfs","psi4i_"+str(part)+"pt"),rhs=BP4.psi4_im_pt[part])],
params="outCverbose=False,CSE_sorting=none") # Generating the CSE for psi4 is the slowest
# operation in this notebook, and much of the CSE
# time is spent sorting CSE expressions. Disabling
# this sorting makes the C codegen 3-4x faster,
# but the tradeoff is that every time this is
# run, the CSE patterns will be different
# (though they should result in mathematically
# *identical* expressions). You can expect
# roundoff-level differences as a result.
with open(os.path.join(Ccodesdir,"Psi4im_pt"+str(part)+"_lowlevel.h"), "w") as file:
file.write(Psi4im_pt)
end = time.time()
print("(BENCH) Finished generating psi4_im_pt"+str(part)+" in "+str(end-start)+" seconds.")
Generating symbolic expressions for psi4... (BENCH) Finished psi4 symbolic expressions in 15.638996601104736 seconds. Output C function psi4() to file BSSN_Two_BHs_Collide_Ccodes_psi4/psi4.h
Instead of measuring $\psi_4$ for all possible (gravitational wave) observers in our simulation domain, we instead decompose it into a natural basis set, which by convention is the spin-weight -2 spherical harmonics.
Here we implement the algorithm for decomposing $\psi_4$ into spin-weight -2 spherical harmonic modes. The decomposition is defined as follows:
$$ {}^{-2}\left[\psi_4\right]_{\ell,m}(t,R) = \int \int \psi_4(t,R,\theta,\phi)\ \left[{}^{-2}Y^*_{\ell,m}(\theta,\phi)\right] \sin \theta d\theta d\phi, $$where
l_max
(set to 2 here) [Back to top]¶Here we output all spin-weight $-2$ spherical harmonics Tutorial Module for $\ell=0$ up to and including $\ell=\ell_{\rm max}$=l_max
(set to 2 here).
import SpinWeight_minus2_SphHarmonics.SpinWeight_minus2_SphHarmonics as swm2
cmd.mkdir(os.path.join(Ccodesdir,"SpinWeight_minus2_SphHarmonics"))
swm2.SpinWeight_minus2_SphHarmonics(maximum_l=l_max,
filename=os.path.join(Ccodesdir,"SpinWeight_minus2_SphHarmonics/SpinWeight_minus2_SphHarmonics.h"))
Note that this diagnostic implementation assumes that Spherical
-like coordinates are used (e.g., SinhSpherical
or Spherical
), which are the most natural coordinate system for decomposing $\psi_4$ into spin-weight -2 modes.
First we process the inputs needed to compute $\psi_4$ at all needed $\theta,\phi$ points
%%writefile $Ccodesdir/driver_psi4_spinweightm2_decomposition.h
void driver_psi4_spinweightm2_decomposition(const paramstruct *restrict params,
const REAL curr_time,const int R_ext_idx,
REAL *restrict xx[3],
const REAL *restrict y_n_gfs,
REAL *restrict diagnostic_output_gfs) {
#include "set_Cparameters.h"
// Step 1: Set the extraction radius R_ext based on the radial index R_ext_idx
REAL R_ext;
{
REAL xx0 = xx[0][R_ext_idx];
REAL xx1 = xx[1][1];
REAL xx2 = xx[2][1];
REAL xCart[3];
xx_to_Cart(params,xx,R_ext_idx,1,1,xCart);
R_ext = sqrt(xCart[0]*xCart[0] + xCart[1]*xCart[1] + xCart[2]*xCart[2]);
}
// Step 2: Compute psi_4 at this extraction radius and store to a local 2D array.
const int sizeof_2Darray = sizeof(REAL)*(Nxx_plus_2NGHOSTS1-2*NGHOSTS)*(Nxx_plus_2NGHOSTS2-2*NGHOSTS);
REAL *restrict psi4r_at_R_ext = (REAL *restrict)malloc(sizeof_2Darray);
REAL *restrict psi4i_at_R_ext = (REAL *restrict)malloc(sizeof_2Darray);
// ... also store theta, sin(theta), and phi to corresponding 1D arrays.
REAL *restrict sinth_array = (REAL *restrict)malloc(sizeof(REAL)*(Nxx_plus_2NGHOSTS1-2*NGHOSTS));
REAL *restrict th_array = (REAL *restrict)malloc(sizeof(REAL)*(Nxx_plus_2NGHOSTS1-2*NGHOSTS));
REAL *restrict ph_array = (REAL *restrict)malloc(sizeof(REAL)*(Nxx_plus_2NGHOSTS2-2*NGHOSTS));
const int i0=R_ext_idx;
#pragma omp parallel for
for(int i1=NGHOSTS;i1<Nxx_plus_2NGHOSTS1-NGHOSTS;i1++) {
th_array[i1-NGHOSTS] = xx[1][i1];
sinth_array[i1-NGHOSTS] = sin(xx[1][i1]);
for(int i2=NGHOSTS;i2<Nxx_plus_2NGHOSTS2-NGHOSTS;i2++) {
ph_array[i2-NGHOSTS] = xx[2][i2];
// Compute real & imaginary parts of psi_4, output to diagnostic_output_gfs
psi4(params, i0,i1,i2, xx, y_n_gfs, diagnostic_output_gfs);
const int idx3d = IDX3S(i0,i1,i2);
const REAL psi4r = (+diagnostic_output_gfs[IDX4ptS(PSI4R_0PTGF,idx3d)]
+diagnostic_output_gfs[IDX4ptS(PSI4R_1PTGF,idx3d)]
+diagnostic_output_gfs[IDX4ptS(PSI4R_2PTGF,idx3d)]);
const REAL psi4i = (+diagnostic_output_gfs[IDX4ptS(PSI4I_0PTGF,idx3d)]
+diagnostic_output_gfs[IDX4ptS(PSI4I_1PTGF,idx3d)]
+diagnostic_output_gfs[IDX4ptS(PSI4I_2PTGF,idx3d)]);
// Store result to "2D" array (actually 1D array with 2D storage):
const int idx2d = (i1-NGHOSTS)*(Nxx_plus_2NGHOSTS2-2*NGHOSTS)+(i2-NGHOSTS);
psi4r_at_R_ext[idx2d] = psi4r;
psi4i_at_R_ext[idx2d] = psi4i;
}
}
Writing BSSN_Two_BHs_Collide_Ccodes_psi4/driver_psi4_spinweightm2_decomposition.h
Next we implement the integral:
$$ {}^{-2}\left[\psi_4\right]_{\ell,m}(t,R) = \int \int \psi_4(t,R,\theta,\phi)\ \left[{}^{-2}Y^*_{\ell,m}(\theta,\phi)\right] \sin \theta d\theta d\phi. $$Since $\psi_4(t,R,\theta,\phi)$ and $\left[{}^{-2}Y^*_{\ell,m}(\theta,\phi)\right]$ are generally complex, for simplicity let's define \begin{align} \psi_4(t,R,\theta,\phi)&=a+i b \\ \left[{}^{-2}Y_{\ell,m}(\theta,\phi)\right] &= c + id\\ \implies \left[{}^{-2}Y^*_{\ell,m}(\theta,\phi)\right] = \left[{}^{-2}Y_{\ell,m}(\theta,\phi)\right]^* &=c-i d \end{align}
Then the product (appearing within the integral) will be given by
\begin{align}
(a + i b) (c-i d) &= (ac + bd) + i(bc - ad),
\end{align}
which cleanly splits the real and complex parts. For better modularity, we output this algorithm to a function decompose_psi4_into_swm2_modes()
in file decompose_psi4_into_swm2_modes.h
. Here, we will call this function from within output_psi4_spinweight_m2_decomposition()
, but in general it could be called from codes that do not use spherical coordinates, and the psi4r_at_R_ext[]
and psi4i_at_R_ext[]
arrays are filled using interpolations.
%%writefile $Ccodesdir/lowlevel_decompose_psi4_into_swm2_modes.h
void lowlevel_decompose_psi4_into_swm2_modes(const paramstruct *restrict params,
const REAL curr_time, const REAL R_ext,
const REAL *restrict th_array,const REAL *restrict sinth_array,const REAL *restrict ph_array,
const REAL *restrict psi4r_at_R_ext,const REAL *restrict psi4i_at_R_ext) {
#include "set_Cparameters.h"
for(int l=2;l<=L_MAX;l++) { // L_MAX is a global variable, since it must be set in Python (so that SpinWeight_minus2_SphHarmonics() computes enough modes)
for(int m=-l;m<=l;m++) {
// Parallelize the integration loop:
REAL psi4r_l_m = 0.0;
REAL psi4i_l_m = 0.0;
#pragma omp parallel for reduction(+:psi4r_l_m,psi4i_l_m)
for(int i1=0;i1<Nxx_plus_2NGHOSTS1-2*NGHOSTS;i1++) {
const REAL th = th_array[i1];
const REAL sinth = sinth_array[i1];
for(int i2=0;i2<Nxx_plus_2NGHOSTS2-2*NGHOSTS;i2++) {
const REAL ph = ph_array[i2];
// Construct integrand for psi4 spin-weight s=-2,l=2,m=0 spherical harmonic
REAL ReY_sm2_l_m,ImY_sm2_l_m;
SpinWeight_minus2_SphHarmonics(l,m, th,ph, &ReY_sm2_l_m,&ImY_sm2_l_m);
const int idx2d = i1*(Nxx_plus_2NGHOSTS2-2*NGHOSTS)+i2;
const REAL a = psi4r_at_R_ext[idx2d];
const REAL b = psi4i_at_R_ext[idx2d];
const REAL c = ReY_sm2_l_m;
const REAL d = ImY_sm2_l_m;
psi4r_l_m += (a*c + b*d) * dxx2 * sinth*dxx1;
psi4i_l_m += (b*c - a*d) * dxx2 * sinth*dxx1;
}
}
// Step 4: Output the result of the integration to file.
char filename[100];
sprintf(filename,"outpsi4_l%d_m%d-%d-r%.2f.txt",l,m, Nxx0,(double)R_ext);
if(m>=0) sprintf(filename,"outpsi4_l%d_m+%d-%d-r%.2f.txt",l,m, Nxx0,(double)R_ext);
FILE *outpsi4_l_m;
// 0 = n*dt when n=0 is exactly represented in double/long double precision,
// so no worries about the result being ~1e-16 in double/ld precision
if(curr_time==0) outpsi4_l_m = fopen(filename, "w");
else outpsi4_l_m = fopen(filename, "a");
fprintf(outpsi4_l_m,"%e %.15e %.15e\n", (double)(curr_time),
(double)psi4r_l_m,(double)psi4i_l_m);
fclose(outpsi4_l_m);
}
}
}
Writing BSSN_Two_BHs_Collide_Ccodes_psi4/lowlevel_decompose_psi4_into_swm2_modes.h
Finally, we complete the function output_psi4_spinweight_m2_decomposition()
, now calling the above routine and freeing all allocated memory.
%%writefile -a $Ccodesdir/driver_psi4_spinweightm2_decomposition.h
// Step 4: Perform integrations across all l,m modes from l=2 up to and including L_MAX (global variable):
lowlevel_decompose_psi4_into_swm2_modes(params, curr_time,R_ext, th_array,sinth_array, ph_array,
psi4r_at_R_ext,psi4i_at_R_ext);
// Step 5: Free all allocated memory:
free(psi4r_at_R_ext); free(psi4i_at_R_ext);
free(sinth_array); free(th_array); free(ph_array);
}
Appending to BSSN_Two_BHs_Collide_Ccodes_psi4/driver_psi4_spinweightm2_decomposition.h
# Step 0: Import the multiprocessing module.
import multiprocessing
# Step 1: Create a list of functions we wish to evaluate in parallel
funcs = [Psi4re,Psi4re,Psi4re,Psi4im,Psi4im,Psi4im, BrillLindquistID,BSSN_RHSs,Ricci,Hamiltonian,gammadet]
# Step 1.a: Define master function for calling all above functions.
# Note that lambdifying this doesn't work in Python 3
def master_func(idx):
if idx < 3: # Call Psi4re(arg)
funcs[idx](idx)
elif idx < 6: # Call Psi4im(arg-3)
funcs[idx](idx-3)
else: # All non-Psi4 functions:
funcs[idx]()
try:
if os.name == 'nt':
# It's a mess to get working in Windows, so we don't bother. :/
# https://medium.com/@grvsinghal/speed-up-your-python-code-using-multiprocessing-on-windows-and-jupyter-or-ipython-2714b49d6fac
raise Exception("Parallel codegen currently not available in Windows")
# Step 1.b: Import the multiprocessing module.
import multiprocessing
# Step 1.c: Evaluate list of functions in parallel if possible;
# otherwise fallback to serial evaluation:
pool = multiprocessing.Pool()
pool.map(master_func,range(len(funcs)))
except:
# Steps 1.b-1.c, alternate: As fallback, evaluate functions in serial.
# This will happen on Android and Windows systems
for idx in range(len(funcs)):
master_func(idx)
Generating C code for psi4_re_pt0 in SinhSpherical coordinates.Generating C code for psi4_re_pt2 in SinhSpherical coordinates.Generating C code for psi4_im_pt1 in SinhSpherical coordinates.Generating C code for psi4_re_pt1 in SinhSpherical coordinates.Generating C code for psi4_im_pt0 in SinhSpherical coordinates. Generating C code for Ricci tensor in SinhSpherical coordinates. Generating C code for BSSN RHSs in SinhSpherical coordinates.Generating C code for psi4_im_pt2 in SinhSpherical coordinates.Generating optimized C code for Brill-Lindquist initial data. May take a while, depending on CoordSystem.Generating optimized C code for gamma constraint. May take a while, depending on CoordSystem. Generating optimized C code for Hamiltonian constraint. May take a while, depending on CoordSystem. Output C function enforce_detgammahat_constraint() to file BSSN_Two_BHs_Collide_Ccodes_psi4/enforce_detgammahat_constraint.h (BENCH) Finished gamma constraint C codegen in 0.13009023666381836 seconds. (BENCH) Finished generating psi4_im_pt1 in 15.4499351978302 seconds. Output C function Ricci_eval() to file BSSN_Two_BHs_Collide_Ccodes_psi4/Ricci_eval.h (BENCH) Finished Ricci C codegen in 21.140114784240723 seconds. (BENCH) Finished generating psi4_im_pt2 in 22.482378721237183 seconds. (BENCH) Finished generating psi4_re_pt2 in 25.301924467086792 seconds. Output C function rhs_eval() to file BSSN_Two_BHs_Collide_Ccodes_psi4/rhs_eval.h (BENCH) Finished BSSN_RHS C codegen in 28.984851121902466 seconds. (BENCH) Finished generating psi4_re_pt1 in 30.5735924243927 seconds. Output C function Hamiltonian_constraint() to file BSSN_Two_BHs_Collide_Ccodes_psi4/Hamiltonian_constraint.h (BENCH) Finished Hamiltonian C codegen in 33.51035833358765 seconds. (BENCH) Finished generating psi4_im_pt0 in 50.90368390083313 seconds. (BENCH) Finished BL initial data codegen in 54.1322066783905 seconds. (BENCH) Finished generating psi4_re_pt0 in 75.61406993865967 seconds.
free_parameters.h
[Back to top]¶Based on declared NRPy+ Cparameters, first we generate declare_Cparameters_struct.h
, set_Cparameters_default.h
, and set_Cparameters[-SIMD].h
.
Then we output free_parameters.h
, which sets initial data parameters, as well as grid domain & reference metric parameters, applying domain_size
and sinh_width
/SymTP_bScale
(if applicable) as set above
# Step 3.f.i: Set free_parameters.h
with open(os.path.join(Ccodesdir,"free_parameters.h"),"w") as file:
file.write("""
// Set free-parameter values.
// Set free-parameter values for BSSN evolution:
params.eta = 2.0;
// Set free parameters for the (Brill-Lindquist) initial data
params.BH1_posn_x = 0.0; params.BH1_posn_y = 0.0; params.BH1_posn_z =+0.25;
params.BH2_posn_x = 0.0; params.BH2_posn_y = 0.0; params.BH2_posn_z =-0.25;
params.BH1_mass = 0.5; params.BH2_mass = 0.5;\n""")
# Append to $Ccodesdir/free_parameters.h reference metric parameters based on generic
# domain_size,sinh_width,sinhv2_const_dr,SymTP_bScale,
# parameters set above.
rfm.out_default_free_parameters_for_rfm(os.path.join(Ccodesdir,"free_parameters.h"),
domain_size,sinh_width,sinhv2_const_dr,SymTP_bScale)
# Step 3.f.ii: Generate set_Nxx_dxx_invdx_params__and__xx.h:
evil_rfm.set_Nxx_dxx_invdx_params__and__xx_h(Ccodesdir)
# Step 3.f.iii: Generate xx_to_Cart.h, which contains xx_to_Cart() for
# (the mapping from xx->Cartesian) for the chosen
# CoordSystem:
evil_rfm.xx_to_Cart_h("xx_to_Cart","./set_Cparameters.h",os.path.join(Ccodesdir,"xx_to_Cart.h"))
# Step 3.f.iv: Generate declare_Cparameters_struct.h, set_Cparameters_default.h, and set_Cparameters[-SIMD].h
import deprecated_NRPy_param_funcs as evil_par
evil_par.generate_Cparameters_Ccodes(os.path.join(Ccodesdir))
Next apply singular, curvilinear coordinate boundary conditions as documented in the corresponding NRPy+ tutorial notebook
import CurviBoundaryConditions.CurviBoundaryConditions as cbcs
cbcs.Set_up_CurviBoundaryConditions(os.path.join(Ccodesdir,"boundary_conditions/"),Cparamspath=os.path.join("../"))
Wrote to file "BSSN_Two_BHs_Collide_Ccodes_psi4/boundary_conditions/parity_conditions_symbolic_dot_products.h" Evolved parity: ( aDD00:4, aDD01:5, aDD02:6, aDD11:7, aDD12:8, aDD22:9, alpha:0, betU0:1, betU1:2, betU2:3, cf:0, hDD00:4, hDD01:5, hDD02:6, hDD11:7, hDD12:8, hDD22:9, lambdaU0:1, lambdaU1:2, lambdaU2:3, trK:0, vetU0:1, vetU1:2, vetU2:3 ) Auxiliary parity: ( H:0, psi4i_0pt:0, psi4i_1pt:0, psi4i_2pt:0, psi4r_0pt:0, psi4r_1pt:0, psi4r_2pt:0 ) AuxEvol parity: ( RbarDD00:4, RbarDD01:5, RbarDD02:6, RbarDD11:7, RbarDD12:8, RbarDD22:9 ) Wrote to file "BSSN_Two_BHs_Collide_Ccodes_psi4/boundary_conditions/EigenCoord_Cart_to_xx.h"
# Part P0: Define REAL, set the number of ghost cells NGHOSTS (from NRPy+'s FD_CENTDERIVS_ORDER),
# and set the CFL_FACTOR (which can be overwritten at the command line)
with open(os.path.join(Ccodesdir,"BSSN_Playground_REAL__NGHOSTS__CFL_FACTOR.h"), "w") as file:
file.write("""
// Part P0.a: Set the number of ghost cells, from NRPy+'s FD_CENTDERIVS_ORDER
#define NGHOSTS """+str(int(FD_order/2)+1)+"""
// Part P0.b: Set the numerical precision (REAL) to double, ensuring all floating point
// numbers are stored to at least ~16 significant digits
#define REAL """+REAL+"""
// Part P0.c: Set the number of ghost cells, from NRPy+'s FD_CENTDERIVS_ORDER
REAL CFL_FACTOR = """+str(CFL_FACTOR)+"""; // Set the CFL Factor. Can be overwritten at command line.
// Part P0.d: We decompose psi_4 into all spin-weight=-2
// l,m spherical harmonics, starting at l=2,
// going up to and including l_max, set here:
#define L_MAX """+str(l_max)+"""
""")
%%writefile $Ccodesdir/BrillLindquist_Playground.c
// Step P0: Define REAL and NGHOSTS; and declare CFL_FACTOR. This header is generated in NRPy+.
#include "BSSN_Playground_REAL__NGHOSTS__CFL_FACTOR.h"
#include "rfm_files/rfm_struct__declare.h"
#include "declare_Cparameters_struct.h"
// All SIMD intrinsics used in SIMD-enabled C code loops are defined here:
#include "SIMD/SIMD_intrinsics.h"
// Step P1: Import needed header files
#include "stdio.h"
#include "stdlib.h"
#include "math.h"
#include "time.h"
#include "stdint.h" // Needed for Windows GCC 6.x compatibility
#ifndef M_PI
#define M_PI 3.141592653589793238462643383279502884L
#endif
#ifndef M_SQRT1_2
#define M_SQRT1_2 0.707106781186547524400844362104849039L
#endif
#define wavespeed 1.0 // Set CFL-based "wavespeed" to 1.0.
// Step P2: Declare the IDX4S(gf,i,j,k) macro, which enables us to store 4-dimensions of
// data in a 1D array. In this case, consecutive values of "i"
// (all other indices held to a fixed value) are consecutive in memory, where
// consecutive values of "j" (fixing all other indices) are separated by
// Nxx_plus_2NGHOSTS0 elements in memory. Similarly, consecutive values of
// "k" are separated by Nxx_plus_2NGHOSTS0*Nxx_plus_2NGHOSTS1 in memory, etc.
#define IDX4S(g,i,j,k) \
( (i) + Nxx_plus_2NGHOSTS0 * ( (j) + Nxx_plus_2NGHOSTS1 * ( (k) + Nxx_plus_2NGHOSTS2 * (g) ) ) )
#define IDX4ptS(g,idx) ( (idx) + (Nxx_plus_2NGHOSTS0*Nxx_plus_2NGHOSTS1*Nxx_plus_2NGHOSTS2) * (g) )
#define IDX3S(i,j,k) ( (i) + Nxx_plus_2NGHOSTS0 * ( (j) + Nxx_plus_2NGHOSTS1 * ( (k) ) ) )
#define LOOP_REGION(i0min,i0max, i1min,i1max, i2min,i2max) \
for(int i2=i2min;i2<i2max;i2++) for(int i1=i1min;i1<i1max;i1++) for(int i0=i0min;i0<i0max;i0++)
#define LOOP_ALL_GFS_GPS(ii) _Pragma("omp parallel for") \
for(int (ii)=0;(ii)<Nxx_plus_2NGHOSTS_tot*NUM_EVOL_GFS;(ii)++)
// Step P3: Set UUGF and VVGF macros, as well as xx_to_Cart()
#include "boundary_conditions/gridfunction_defines.h"
// Step P4: Set xx_to_Cart(const paramstruct *restrict params,
// REAL *restrict xx[3],
// const int i0,const int i1,const int i2,
// REAL xCart[3]),
// which maps xx->Cartesian via
// {xx[0][i0],xx[1][i1],xx[2][i2]}->{xCart[0],xCart[1],xCart[2]}
#include "xx_to_Cart.h"
// Step P5: Defines set_Nxx_dxx_invdx_params__and__xx(const int EigenCoord, const int Nxx[3],
// paramstruct *restrict params, REAL *restrict xx[3]),
// which sets params Nxx,Nxx_plus_2NGHOSTS,dxx,invdx, and xx[] for
// the chosen Eigen-CoordSystem if EigenCoord==1, or
// CoordSystem if EigenCoord==0.
#include "set_Nxx_dxx_invdx_params__and__xx.h"
// Step P6: Include basic functions needed to impose curvilinear
// parity and boundary conditions.
#include "boundary_conditions/CurviBC_include_Cfunctions.h"
// Step P7: Implement the algorithm for upwinding.
// *NOTE*: This upwinding is backwards from
// usual upwinding algorithms, because the
// upwinding control vector in BSSN (the shift)
// acts like a *negative* velocity.
//#define UPWIND_ALG(UpwindVecU) UpwindVecU > 0.0 ? 1.0 : 0.0
// Step P8: Include function for enforcing detgammabar constraint.
#include "enforce_detgammahat_constraint.h"
// Step P9: Find the CFL-constrained timestep
#include "find_timestep.h"
// Step P10: Declare function necessary for setting up the initial data.
// Step P10.a: Define BSSN_ID() for BrillLindquist initial data
// Step P10.b: Set the generic driver function for setting up BSSN initial data
#include "initial_data.h"
// Step P11: Declare function for evaluating Hamiltonian constraint (diagnostic)
#include "Hamiltonian_constraint.h"
// Step P12: Declare rhs_eval function, which evaluates BSSN RHSs
#include "rhs_eval.h"
// Step P13: Declare Ricci_eval function, which evaluates Ricci tensor
#include "Ricci_eval.h"
// Step P14: Declare function for evaluating real and imaginary parts of psi4 (diagnostic)
#include "psi4.h"
#include "SpinWeight_minus2_SphHarmonics/SpinWeight_minus2_SphHarmonics.h"
#include "lowlevel_decompose_psi4_into_swm2_modes.h"
#include "driver_psi4_spinweightm2_decomposition.h"
// main() function:
// Step 0: Read command-line input, set up grid structure, allocate memory for gridfunctions, set up coordinates
// Step 1: Set up initial data to an exact solution
// Step 2: Start the timer, for keeping track of how fast the simulation is progressing.
// Step 3: Integrate the initial data forward in time using the chosen RK-like Method of
// Lines timestepping algorithm, and output periodic simulation diagnostics
// Step 3.a: Output 2D data file periodically, for visualization
// Step 3.b: Step forward one timestep (t -> t+dt) in time using
// chosen RK-like MoL timestepping algorithm
// Step 3.c: If t=t_final, output conformal factor & Hamiltonian
// constraint violation to 2D data file
// Step 3.d: Progress indicator printing to stderr
// Step 4: Free all allocated memory
int main(int argc, const char *argv[]) {
paramstruct params;
#include "set_Cparameters_default.h"
// Step 0a: Read command-line input, error out if nonconformant
if((argc != 4 && argc != 5) || atoi(argv[1]) < NGHOSTS || atoi(argv[2]) < NGHOSTS || atoi(argv[3]) < 2 /* FIXME; allow for axisymmetric sims */) {
fprintf(stderr,"Error: Expected three command-line arguments: ./BrillLindquist_Playground Nx0 Nx1 Nx2,\n");
fprintf(stderr,"where Nx[0,1,2] is the number of grid points in the 0, 1, and 2 directions.\n");
fprintf(stderr,"Nx[] MUST BE larger than NGHOSTS (= %d)\n",NGHOSTS);
exit(1);
}
if(argc == 5) {
CFL_FACTOR = strtod(argv[4],NULL);
}
if(CFL_FACTOR > 0.5 && atoi(argv[3])!=2) {
fprintf(stderr,"WARNING: CFL_FACTOR was set to %e, which is > 0.5.\n",CFL_FACTOR);
fprintf(stderr," This will generally only be stable if the simulation is purely axisymmetric\n");
fprintf(stderr," However, Nx2 was set to %d>2, which implies a non-axisymmetric simulation\n",atoi(argv[3]));
}
// Step 0b: Set up numerical grid structure, first in space...
const int Nxx[3] = { atoi(argv[1]), atoi(argv[2]), atoi(argv[3]) };
if(Nxx[0]%2 != 0 || Nxx[1]%2 != 0 || Nxx[2]%2 != 0) {
fprintf(stderr,"Error: Cannot guarantee a proper cell-centered grid if number of grid cells not set to even number.\n");
fprintf(stderr," For example, in case of angular directions, proper symmetry zones will not exist.\n");
exit(1);
}
// Step 0c: Set free parameters, overwriting Cparameters defaults
// by hand or with command-line input, as desired.
#include "free_parameters.h"
// Step 0d: Uniform coordinate grids are stored to *xx[3]
REAL *xx[3];
// Step 0d.i: Set bcstruct
bc_struct bcstruct;
{
int EigenCoord = 1;
// Step 0d.ii: Call set_Nxx_dxx_invdx_params__and__xx(), which sets
// params Nxx,Nxx_plus_2NGHOSTS,dxx,invdx, and xx[] for the
// chosen Eigen-CoordSystem.
set_Nxx_dxx_invdx_params__and__xx(EigenCoord, Nxx, ¶ms, xx);
// Step 0d.iii: Set Nxx_plus_2NGHOSTS_tot
#include "set_Cparameters-nopointer.h"
const int Nxx_plus_2NGHOSTS_tot = Nxx_plus_2NGHOSTS0*Nxx_plus_2NGHOSTS1*Nxx_plus_2NGHOSTS2;
// Step 0e: Find ghostzone mappings; set up bcstruct
#include "boundary_conditions/driver_bcstruct.h"
// Step 0e.i: Free allocated space for xx[][] array
for(int i=0;i<3;i++) free(xx[i]);
}
// Step 0f: Call set_Nxx_dxx_invdx_params__and__xx(), which sets
// params Nxx,Nxx_plus_2NGHOSTS,dxx,invdx, and xx[] for the
// chosen (non-Eigen) CoordSystem.
int EigenCoord = 0;
set_Nxx_dxx_invdx_params__and__xx(EigenCoord, Nxx, ¶ms, xx);
// Step 0g: Set all C parameters "blah" for params.blah, including
// Nxx_plus_2NGHOSTS0 = params.Nxx_plus_2NGHOSTS0, etc.
#include "set_Cparameters-nopointer.h"
const int Nxx_plus_2NGHOSTS_tot = Nxx_plus_2NGHOSTS0*Nxx_plus_2NGHOSTS1*Nxx_plus_2NGHOSTS2;
// Step 0h: Time coordinate parameters
const REAL t_final = 1.5*domain_size; /* Final time is set so that at t=t_final,
* data at the origin have not been corrupted
* by the approximate outer boundary condition */
// Step 0i: Set timestep based on smallest proper distance between gridpoints and CFL factor
REAL dt = find_timestep(¶ms, xx);
//fprintf(stderr,"# Timestep set to = %e\n",(double)dt);
int N_final = (int)(t_final / dt + 0.5); // The number of points in time.
// Add 0.5 to account for C rounding down
// typecasts to integers.
const REAL out_approx_every_t = 0.5;
int output_every_N = (int)(out_approx_every_t*((REAL)N_final)/t_final);
// Step 0j: Error out if the number of auxiliary gridfunctions outnumber evolved gridfunctions.
// This is a limitation of the RK method. You are always welcome to declare & allocate
// additional gridfunctions by hand.
if(NUM_AUX_GFS > NUM_EVOL_GFS) {
fprintf(stderr,"Error: NUM_AUX_GFS > NUM_EVOL_GFS. Either reduce the number of auxiliary gridfunctions,\n");
fprintf(stderr," or allocate (malloc) by hand storage for *diagnostic_output_gfs. \n");
exit(1);
}
// Step 0k: Allocate memory for gridfunctions
#include "MoLtimestepping/RK_Allocate_Memory.h"
REAL *restrict auxevol_gfs = (REAL *)malloc(sizeof(REAL) * NUM_AUXEVOL_GFS * Nxx_plus_2NGHOSTS_tot);
// Step 0l: Set up precomputed reference metric arrays
// Step 0l.i: Allocate space for precomputed reference metric arrays.
#include "rfm_files/rfm_struct__malloc.h"
// Step 0l.ii: Define precomputed reference metric arrays.
{
#include "set_Cparameters-nopointer.h"
#include "rfm_files/rfm_struct__define.h"
}
// Step 1: Set up initial data to an exact solution
initial_data(¶ms, xx, y_n_gfs);
// Step 1b: Apply boundary conditions, as initial data
// are sometimes ill-defined in ghost zones.
// E.g., spherical initial data might not be
// properly defined at points where r=-1.
apply_bcs_curvilinear(¶ms, &bcstruct, NUM_EVOL_GFS,evol_gf_parity, y_n_gfs);
enforce_detgammahat_constraint(&rfmstruct, ¶ms, y_n_gfs);
// Step 2: Start the timer, for keeping track of how fast the simulation is progressing.
#ifdef __linux__ // Use high-precision timer in Linux.
struct timespec start, end;
clock_gettime(CLOCK_REALTIME, &start);
#else // Resort to low-resolution, standards-compliant timer in non-Linux OSs
// http://www.cplusplus.com/reference/ctime/time/
time_t start_timer,end_timer;
time(&start_timer); // Resolution of one second...
#endif
// Step 3: Integrate the initial data forward in time using the chosen RK-like Method of
// Lines timestepping algorithm, and output periodic simulation diagnostics
for(int n=0;n<=N_final;n++) { // Main loop to progress forward in time.
// Step 3.a: Output diagnostics
if(n%output_every_N == 0) {
// Step 3.a.i: Output psi4 spin-weight -2 decomposed data, every N_output_every
for(int r_ext_idx = (Nxx_plus_2NGHOSTS0-NGHOSTS)/4;
r_ext_idx<(Nxx_plus_2NGHOSTS0-NGHOSTS)*0.9;
r_ext_idx+=5) {
// psi_4 mode-by-mode spin-weight -2 spherical harmonic decomposition routine
driver_psi4_spinweightm2_decomposition(¶ms, ((REAL)n)*dt,r_ext_idx,
xx, y_n_gfs, diagnostic_output_gfs);
}
}
if(n%100 == 0) {
// Step 3.a.ii: Regularly output conformal factor & Hamiltonian
// constraint violation to 2D data file
// Evaluate Hamiltonian constraint violation
Hamiltonian_constraint(&rfmstruct, ¶ms, y_n_gfs, diagnostic_output_gfs);
char filename[100];
sprintf(filename,"out%d-%08d.txt",Nxx[0],n);
FILE *out2D = fopen(filename, "w");
LOOP_REGION(NGHOSTS,Nxx_plus_2NGHOSTS0-NGHOSTS,
NGHOSTS,Nxx_plus_2NGHOSTS1-NGHOSTS,
NGHOSTS,Nxx_plus_2NGHOSTS2-NGHOSTS) {
REAL xCart[3]; xx_to_Cart(¶ms,xx,i0,i1,i2, xCart);
int idx = IDX3S(i0,i1,i2);
fprintf(out2D,"%e %e %e %e\n",xCart[1],xCart[2],
y_n_gfs[IDX4ptS(CFGF,idx)], log10(fabs(diagnostic_output_gfs[IDX4ptS(HGF,idx)])));
}
fclose(out2D);
}
// Step 3.b: Step forward one timestep (t -> t+dt) in time using
// chosen RK-like MoL timestepping algorithm
#include "MoLtimestepping/RK_MoL.h"
// Step 3.c: If t=t_final, output conformal factor & Hamiltonian
// constraint violation to 2D data file
if(n==N_final-1) {
// Evaluate Hamiltonian constraint violation
Hamiltonian_constraint(&rfmstruct, ¶ms, y_n_gfs, diagnostic_output_gfs);
char filename[100];
sprintf(filename,"out%d.txt",Nxx[0]);
FILE *out2D = fopen(filename, "w");
LOOP_REGION(NGHOSTS,Nxx_plus_2NGHOSTS0-NGHOSTS,
NGHOSTS,Nxx_plus_2NGHOSTS1-NGHOSTS,
NGHOSTS,Nxx_plus_2NGHOSTS2-NGHOSTS) {
REAL xCart[3]; xx_to_Cart(¶ms,xx,i0,i1,i2, xCart);
int idx = IDX3S(i0,i1,i2);
fprintf(out2D,"%e %e %e %e\n",xCart[1],xCart[2], y_n_gfs[IDX4ptS(CFGF,idx)],
log10(fabs(diagnostic_output_gfs[IDX4ptS(HGF,idx)])));
}
fclose(out2D);
}
// Step 3.d: Progress indicator printing to stderr
// Step 3.d.i: Measure average time per iteration
#ifdef __linux__ // Use high-precision timer in Linux.
clock_gettime(CLOCK_REALTIME, &end);
const long long unsigned int time_in_ns = 1000000000L * (end.tv_sec - start.tv_sec) + end.tv_nsec - start.tv_nsec;
#else // Resort to low-resolution, standards-compliant timer in non-Linux OSs
time(&end_timer); // Resolution of one second...
REAL time_in_ns = difftime(end_timer,start_timer)*1.0e9+0.5; // Round up to avoid divide-by-zero.
#endif
const REAL s_per_iteration_avg = ((REAL)time_in_ns / (REAL)n) / 1.0e9;
const int iterations_remaining = N_final - n;
const REAL time_remaining_in_mins = s_per_iteration_avg * (REAL)iterations_remaining / 60.0;
const REAL num_RHS_pt_evals = (REAL)(Nxx[0]*Nxx[1]*Nxx[2]) * 4.0 * (REAL)n; // 4 RHS evals per gridpoint for RK4
const REAL RHS_pt_evals_per_sec = num_RHS_pt_evals / ((REAL)time_in_ns / 1.0e9);
// Step 3.d.ii: Output simulation progress to stderr
if(n % 40 == 0) {
fprintf(stderr,"%c[2K", 27); // Clear the line
fprintf(stderr,"It: %d t=%.2f dt=%.2e | %.1f%%; ETA %.0f s | t/h %.2f | gp/s %.2e\r", // \r is carriage return, move cursor to the beginning of the line
n, n * (double)dt, (double)dt, (double)(100.0 * (REAL)n / (REAL)N_final),
(double)time_remaining_in_mins*60, (double)(dt * 3600.0 / s_per_iteration_avg), (double)RHS_pt_evals_per_sec);
fflush(stderr); // Flush the stderr buffer
} // End progress indicator if(n % 10 == 0)
} // End main loop to progress forward in time.
fprintf(stderr,"\n"); // Clear the final line of output from progress indicator.
// Step 4: Free all allocated memory
#include "rfm_files/rfm_struct__freemem.h"
#include "boundary_conditions/bcstruct_freemem.h"
#include "MoLtimestepping/RK_Free_Memory.h"
free(auxevol_gfs);
for(int i=0;i<3;i++) free(xx[i]);
return 0;
}
Writing BSSN_Two_BHs_Collide_Ccodes_psi4/BrillLindquist_Playground.c
if EvolOption == "low resolution":
Nr = 270
Ntheta = 8
elif EvolOption == "high resolution":
Nr = 800
Ntheta = 16
else:
print("Error: unknown EvolOption!")
sys.exit(1)
import cmdline_helper as cmd
print("Now compiling, should take ~10 seconds...\n")
start = time.time()
cmd.C_compile(os.path.join(Ccodesdir, "BrillLindquist_Playground.c"),
os.path.join(outdir, "BrillLindquist_Playground"),
compile_mode="optimized")
end = time.time()
print("(BENCH) Finished in "+str(end-start)+" seconds.\n\n")
if Nr == 800:
print("Now running. Should take ~8 hours...\n")
if Nr == 270:
print("Now running. Should take ~30 minutes...\n")
# Change to output directory
os.chdir(os.path.join(outdir))
cmd.delete_existing_files("out*.txt")
cmd.delete_existing_files("out*.png")
cmd.Execute("BrillLindquist_Playground", str(Nr)+" "+str(Ntheta)+" 2")
# Change back to root directory
os.chdir(os.path.join("..", ".."))
Now compiling, should take ~10 seconds... Compiling executable... (EXEC): Executing `gcc -std=gnu99 -Ofast -fopenmp -march=native -funroll-loops BSSN_Two_BHs_Collide_Ccodes_psi4/BrillLindquist_Playground.c -o BSSN_Two_BHs_Collide_Ccodes_psi4/output/BrillLindquist_Playground -lm`... (BENCH): Finished executing in 28.450400590896606 seconds. Finished compilation. (BENCH) Finished in 28.45347309112549 seconds. Now running. Should take ~30 minutes... (EXEC): Executing `taskset -c 1,3,5,7,9,11,13,15 ./BrillLindquist_Playground 270 8 2`... It: 30600 t=224.92 dt=7.35e-03 | 100.0%; ETA 0 s | t/h 1435.93 | gp/s 9.38e+055 (BENCH): Finished executing in 564.2851254940033 seconds.
According to black hole perturbation theory (Berti et al), the resultant black hole should ring down with dominant, spin-weight $s=-2$ spherical harmonic mode $(l=2,m=0)$ according to
$$ {}_{s=-2}\text{Re}(\psi_4)_{l=2,m=0} = A e^{−0.0890 t/M} \cos(0.3737 t/M+ \phi), $$where $M=1$ for these data, and $A$ and $\phi$ are an arbitrary amplitude and phase, respectively. Here we will plot the resulting waveform at $r/M=33.13$, comparing to the expected frequency and amplitude falloff predicted by black hole perturbation theory.
Notice that we find about 4.2 orders of magnitude agreement! If you are willing to invest more resources and wait much longer, you will find approximately 8.5 orders of magnitude agreement (better than Fig 6 of Ruchlin et al) if you adjust the above code parameters such that
AMPL
) set to 300t_final
) set to 275BH1_posn_z = -BH2_posn_z = 0.25
%matplotlib inline
import numpy as np
import matplotlib.pyplot as plt
from matplotlib.pyplot import savefig
from matplotlib import rc
rc('text', usetex=True)
if Nr == 270:
extraction_radius = "30.20"
Amplitude = 1.8e-2
Phase = 2.8
elif Nr == 800:
extraction_radius = "33.64"
Amplitude = 1.8e-2
Phase = 2.8
else:
print("Error: output is not tuned for Nr = "+str(Nr)+" . Plotting disabled.")
exit(1)
# Transposed for easier unpacking:
t,psi4r,psi4i = np.loadtxt(os.path.join(outdir, "outpsi4_l2_m+0-"+str(Nr)+"-r"+extraction_radius+".txt")).T
t_retarded = []
log10abspsi4r = []
bh_pert_thry = []
for i in range(len(psi4r)):
retarded_time = t[i]-float(extraction_radius)
t_retarded.append(retarded_time)
log10abspsi4r.append(np.log(float(extraction_radius)*np.abs(psi4r[i]))/np.log(10))
bh_pert_thry.append(np.log(Amplitude*np.exp(-0.0890*retarded_time)*np.abs(np.cos(0.3737*retarded_time+Phase)))/np.log(10))
# print(bh_pert_thry)
fig, ax = plt.subplots()
plt.title("Grav. Wave Agreement with BH perturbation theory",fontsize=18)
plt.xlabel("$(t - R_{ext})/M$",fontsize=16)
plt.ylabel('$\log_{10}|\psi_4|$',fontsize=16)
ax.plot(t_retarded, log10abspsi4r, 'k-', label='SENR/NRPy+ simulation')
ax.plot(t_retarded, bh_pert_thry, 'k--', label='BH perturbation theory')
# ax.set_xlim([0,t_retarded[len(psi4r1)-1]])
ax.set_xlim([0, 1.5*domain_size - float(extraction_radius)])
ax.set_ylim([-13.5, -1.5])
plt.xticks(size = 14)
plt.yticks(size = 14)
legend = ax.legend(loc='upper right', shadow=True, fontsize='x-large')
legend.get_frame().set_facecolor('C1')
plt.show()
# Note that you'll need `dvipng` installed to generate the following file:
savefig(os.path.join(outdir, "BHperttheorycompare.png"), dpi=150)
<Figure size 432x288 with 0 Axes>
The following code cell converts this Jupyter notebook into a proper, clickable $\LaTeX$-formatted PDF file. After the cell is successfully run, the generated PDF may be found in the root NRPy+ tutorial directory, with filename Tutorial-Start_to_Finish-BSSNCurvilinear-Two_BHs_Collide-Psi4.pdf (Note that clicking on this link may not work; you may need to open the PDF file through another means.)
import cmdline_helper as cmd # NRPy+: Multi-platform Python command-line interface
cmd.output_Jupyter_notebook_to_LaTeXed_PDF("Tutorial-Start_to_Finish-BSSNCurvilinear-Two_BHs_Collide-Psi4")
Created Tutorial-Start_to_Finish-BSSNCurvilinear-Two_BHs_Collide-Psi4.tex, and compiled LaTeX file to PDF file Tutorial-Start_to_Finish- BSSNCurvilinear-Two_BHs_Collide-Psi4.pdf