def BrillLindquist(ComputeADMGlobalsOnly=False, include_NRPy_basic_defines_and_pickle=False): # Step 2: Setting up Brill-Lindquist initial data # Step 2.a: Set spatial dimension (must be 3 for BSSN) DIM = 3 par.set_parval_from_str("grid::DIM", DIM) global Cartxyz, gammaCartDD, KCartDD, alphaCart, betaCartU, BCartU Cartxyz = ixp.declarerank1("Cartxyz") # Step 2.b: Set psi, the conformal factor: psi = sp.sympify(1) psi += BH1_mass / (2 * sp.sqrt((Cartxyz[0] - BH1_posn_x)**2 + (Cartxyz[1] - BH1_posn_y)**2 + (Cartxyz[2] - BH1_posn_z)**2)) psi += BH2_mass / (2 * sp.sqrt((Cartxyz[0] - BH2_posn_x)**2 + (Cartxyz[1] - BH2_posn_y)**2 + (Cartxyz[2] - BH2_posn_z)**2)) # Step 2.c: Set all needed ADM variables in Cartesian coordinates gammaCartDD = ixp.zerorank2() KCartDD = ixp.zerorank2() # K_{ij} = 0 for these initial data for i in range(DIM): gammaCartDD[i][i] = psi**4 alphaCart = 1 / psi**2 betaCartU = ixp.zerorank1( ) # We generally choose \beta^i = 0 for these initial data BCartU = ixp.zerorank1( ) # We generally choose B^i = 0 for these initial data if ComputeADMGlobalsOnly == True: return cf,hDD,lambdaU,aDD,trK,alpha,vetU,betU = \ AtoB.Convert_Spherical_or_Cartesian_ADM_to_BSSN_curvilinear("Cartesian",Cartxyz, gammaCartDD,KCartDD,alphaCart,betaCartU,BCartU) import BSSN.BSSN_ID_function_string as bIDf # Generates initial_data() C function & stores to outC_function_dict["initial_data"] bIDf.BSSN_ID_function_string( cf, hDD, lambdaU, aDD, trK, alpha, vetU, betU, include_NRPy_basic_defines=include_NRPy_basic_defines_and_pickle) if include_NRPy_basic_defines_and_pickle: return pickle_NRPy_env()
def UIUCBlackHole(ComputeADMGlobalsOnly=False, include_NRPy_basic_defines_and_pickle=False): global Sph_r_th_ph, r, th, ph, gammaSphDD, KSphDD, alphaSph, betaSphU, BSphU # All gridfunctions will be written in terms of spherical coordinates (r, th, ph): r, th, ph = sp.symbols('r th ph', real=True) # Step 0: Set spatial dimension (must be 3 for BSSN) DIM = 3 par.set_parval_from_str("grid::DIM", DIM) # Step 1: Set psi, the conformal factor: # Spin per unit mass a = M * chi # Defined under equation 1 in Liu, Etienne, & Shapiro (2009) https://arxiv.org/pdf/1001.4077.pdf # Boyer - Lindquist outer horizon rp = M + sp.sqrt(M**2 - a**2) # Boyer - Lindquist inner horizon rm = M - sp.sqrt(M**2 - a**2) # Boyer - Lindquist radius in terms of UIUC radius # Eq. 11 # r_{BL} = r * ( 1 + r_+ / 4r )^2 rBL = r * (1 + rp / (4 * r))**2 # Expressions found below Eq. 2 # Sigma = r_{BL}^2 + a^2 cos^2 theta SIG = rBL**2 + a**2 * sp.cos(th)**2 # Delta = r_{BL}^2 - 2Mr_{BL} + a^2 DEL = rBL**2 - 2 * M * rBL + a**2 # A = (r_{BL}^2 + a^2)^2 - Delta a^2 sin^2 theta AA = (rBL**2 + a**2)**2 - DEL * a**2 * sp.sin(th)**2 # *** The ADM 3-metric in spherical basis *** gammaSphDD = ixp.zerorank2() # Declare the nonzero components of the 3-metric # (Eq. 13 of Liu, Etienne, & Shapiro, https://arxiv.org/pdf/1001.4077.pdf): # ds^2 = Sigma (r + r_+/4)^2 / ( r^3 (r_{BL} - r_- ) * dr^2 + # Sigma d theta^2 + (A sin^2 theta) / Sigma * d\phi^2 gammaSphDD[0][0] = ((SIG * (r + rp / 4)**2) / (r**3 * (rBL - rm))) gammaSphDD[1][1] = SIG gammaSphDD[2][2] = AA / SIG * sp.sin(th)**2 # *** The physical trace-free extrinsic curvature in spherical basis *** # Nonzero components of the extrinsic curvature K, given by # Eq. 14 of Liu, Etienne, & Shapiro, https://arxiv.org/pdf/1001.4077.pdf: KSphDD = ixp.zerorank2() # K_{r phi} = K_{phi r} = (Ma sin^2 theta) / (Sigma sqrt{A Sigma}) * # [3r^4_{BL} + 2a^2 r^2_{BL} - a^4 - a^2 (r^2_{BL} - a^2) sin^2 theta] * # (1 + r_+ / 4r) (1 / sqrt{r(r_{BL} - r_-)}) KSphDD[0][2] = KSphDD[2][0] = (M*a*sp.sin(th)**2)/(SIG*sp.sqrt(AA*SIG))*\ (3*rBL**4 + 2*a**2*rBL**2 - a**4- a**2*(rBL**2 - a**2)*\ sp.sin(th)**2)*(1 + rp/(4*r))*1/sp.sqrt(r*(rBL - rm)) # Components of the extrinsic curvature K, given by # Eq. 15 of Liu, Etienne, & Shapiro, https://arxiv.org/pdf/1001.4077.pdf: # K_{theta phi} = K_{phi theta} = -(2a^3 Mr_{BL} cos theta sin^3 theta) / # (Sigma sqrt{A Sigma}) x (r - r_+ / 4) sqrt{(r_{BL} - r_-) / r } KSphDD[1][2] = KSphDD[2][1] = -((2*a**3*M*rBL*sp.cos(th)*sp.sin(th)**3)/ \ (SIG*sp.sqrt(AA*SIG)))*(r - rp/4)*sp.sqrt((rBL - rm)/r) alphaSph = sp.sympify( 1 ) # We generally choose alpha = 1/psi**2 (psi = BSSN conformal factor) for these initial data betaSphU = ixp.zerorank1( ) # We generally choose \beta^i = 0 for these initial data BSphU = ixp.zerorank1( ) # We generally choose B^i = 0 for these initial data if ComputeADMGlobalsOnly: return # Validated against original SENR: KSphDD[0][2], KSphDD[1][2], gammaSphDD[2][2], gammaSphDD[0][0], gammaSphDD[1][1] #print(sp.mathematica_code(gammaSphDD[1][1])) Sph_r_th_ph = [r, th, ph] cf,hDD,lambdaU,aDD,trK,alpha,vetU,betU = \ AtoB.Convert_Spherical_or_Cartesian_ADM_to_BSSN_curvilinear("Spherical", Sph_r_th_ph, gammaSphDD,KSphDD,alphaSph,betaSphU,BSphU) # Let's choose alpha = 1/psi**2 (psi = BSSN conformal factor) for these initial data, # where psi = exp(phi); chi = 1/psi**4; W = 1/psi**2 if par.parval_from_str("EvolvedConformalFactor_cf") == "phi": alpha = sp.exp(-2 * cf) elif par.parval_from_str("EvolvedConformalFactor_cf") == "chi": alpha = sp.sqrt(cf) elif par.parval_from_str("EvolvedConformalFactor_cf") == "W": alpha = cf else: print("Error EvolvedConformalFactor_cf type = \"" + par.parval_from_str("EvolvedConformalFactor_cf") + "\" unknown.") sys.exit(1) import BSSN.BSSN_ID_function_string as bIDf # Generates initial_data() C function & stores to outC_function_dict["initial_data"] bIDf.BSSN_ID_function_string( cf, hDD, lambdaU, aDD, trK, alpha, vetU, betU, include_NRPy_basic_defines=include_NRPy_basic_defines_and_pickle) if include_NRPy_basic_defines_and_pickle: return pickle_NRPy_env()
def add_psi4_tetrad_to_Cfunction_dict(includes=None, rel_path_to_Cparams=os.path.join("."), setPsi4tozero=False): starttime = print_msg_with_timing("psi4 tetrads", msg="Ccodegen", startstop="start") # Set up the C function for BSSN basis transformations desc = "Compute tetrad for psi4" name = "psi4_tetrad" # First set up the symbolic expressions (RHSs) and their names (LHSs) psi4tet.Psi4_tetrads() list_of_varnames = [] list_of_symbvars = [] for i in range(4): list_of_varnames.append("*mre4U" + str(i)) list_of_symbvars.append(psi4tet.mre4U[i]) for i in range(4): list_of_varnames.append("*mim4U" + str(i)) list_of_symbvars.append(psi4tet.mim4U[i]) for i in range(4): list_of_varnames.append("*n4U" + str(i)) list_of_symbvars.append(psi4tet.n4U[i]) paramsindent = " " params = """const paramstruct *restrict params,\n""" + paramsindent list_of_metricvarnames = ["cf"] for i in range(3): for j in range(i, 3): list_of_metricvarnames.append("hDD" + str(i) + str(j)) for var in list_of_metricvarnames: params += "const REAL " + var + "," params += "\n" + paramsindent for var in list_of_varnames: params += "REAL " + var + "," params += "\n" + paramsindent + "REAL *restrict xx[3], const int i0,const int i1,const int i2" # Set the body of the function body = "" outCparams = "includebraces=False,outCverbose=False,enable_SIMD=False,preindent=1" if not setPsi4tozero: for i in range(3): body += " const REAL xx" + str(i) + " = xx[" + str( i) + "][i" + str(i) + "];\n" body += " // Compute tetrads:\n" body += " {\n" # Sort the lhss list alphabetically, and rhss to match: lhss, rhss = [ list(x) for x in zip(*sorted(zip(list_of_varnames, list_of_symbvars), key=lambda pair: pair[0])) ] body += outputC(rhss, lhss, filename="returnstring", params=outCparams) body += " }\n" elif setPsi4tozero: body += "return;\n" loopopts = "" print_msg_with_timing("psi4 tetrads", msg="Ccodegen", startstop="stop", starttime=starttime) add_to_Cfunction_dict(includes=includes, desc=desc, name=name, params=params, body=body, loopopts=loopopts, rel_path_to_Cparams=rel_path_to_Cparams) return pickle_NRPy_env()
def add_SpinWeight_minus2_SphHarmonics_to_Cfunction_dict( includes=None, rel_path_to_Cparams=os.path.join("."), maximum_l=8): starttime = print_msg_with_timing("Spin-weight s=-2 Spherical Harmonics", msg="Ccodegen", startstop="start") # Set up the C function for computing the spin-weight -2 spherical harmonic at theta,phi: Y_{s=-2, l,m}(theta,phi) prefunc = r"""// Compute at a single point (th,ph) the spin-weight -2 spherical harmonic Y_{s=-2, l,m}(th,ph) // Manual "inline void" of this function results in compilation error with clang. void SpinWeight_minus2_SphHarmonics(const int l, const int m, const REAL th, const REAL ph, REAL *reYlmswm2_l_m, REAL *imYlmswm2_l_m) { """ # Construct prefunc: outCparams = "preindent=1,outCfileaccess=a,outCverbose=False,includebraces=False" prefunc += """ switch(l) { """ for l in range(maximum_l + 1): # Output values up to and including l=8. prefunc += " case " + str(l) + ":\n" prefunc += " switch(m) {\n" for m in range(-l, l + 1): prefunc += " case " + str(m) + ":\n" prefunc += " {\n" Y_m2_lm = SWm2SH.Y(-2, l, m, SWm2SH.th, SWm2SH.ph) prefunc += outputC([sp.re(Y_m2_lm), sp.im(Y_m2_lm)], ["*reYlmswm2_l_m", "*imYlmswm2_l_m"], "returnstring", outCparams) prefunc += " }\n" prefunc += " return;\n" prefunc += " } // END switch(m)\n" prefunc += " } // END switch(l)\n" prefunc += r""" fprintf(stderr, "ERROR: SpinWeight_minus2_SphHarmonics handles only l=[0,""" + str( maximum_l) + r"""] and only m=[-l,+l] is defined.\n"); fprintf(stderr, " You chose l=%d and m=%d, which is out of these bounds.\n",l,m); exit(1); } void lowlevel_decompose_psi4_into_swm2_modes(const int Nxx_plus_2NGHOSTS1,const int Nxx_plus_2NGHOSTS2, const REAL dxx1, const REAL dxx2, 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) { for(int l=2;l<=""" + str( maximum_l) + r""";l++) { // The maximum l here is set in Python. 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-r%.2f.txt",l,m, (double)R_ext); // If you love "+"'s in filenames by all means enable this (ugh): //if(m>=0) sprintf(filename,"outpsi4_l%d_m+%d-r%.2f.txt",l,m, (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); } } } """ desc = "" name = "driver__spherlikegrids__psi4_spinweightm2_decomposition" params = r"""const paramstruct *restrict params, REAL *restrict diagnostic_output_gfs, const int *restrict list_of_R_ext_idxs, const int num_of_R_ext_idxs, const REAL time, REAL *restrict xx[3],void xx_to_Cart(const paramstruct *restrict params, REAL *restrict xx[3],const int i0,const int i1,const int i2, REAL xCart[3])""" body = r""" // Step 1: Allocate memory for 2D arrays used to store psi4, theta, sin(theta), and phi. 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)); // Step 2: Loop over all extraction indices: for(int ii0=0;ii0<num_of_R_ext_idxs;ii0++) { // Step 2.a: Set the extraction radius R_ext based on the radial index R_ext_idx REAL R_ext; { REAL xCart[3]; xx_to_Cart(params,xx,list_of_R_ext_idxs[ii0],1,1,xCart); // values for itheta and iphi don't matter. R_ext = sqrt(xCart[0]*xCart[0] + xCart[1]*xCart[1] + xCart[2]*xCart[2]); } // Step 2.b: Compute psi_4 at this extraction radius and store to a local 2D array. const int i0=list_of_R_ext_idxs[ii0]; #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 const REAL psi4r = (diagnostic_output_gfs[IDX4S(PSI4_PART0REGF, i0,i1,i2)] + diagnostic_output_gfs[IDX4S(PSI4_PART1REGF, i0,i1,i2)] + diagnostic_output_gfs[IDX4S(PSI4_PART2REGF, i0,i1,i2)]); const REAL psi4i = (diagnostic_output_gfs[IDX4S(PSI4_PART0IMGF, i0,i1,i2)] + diagnostic_output_gfs[IDX4S(PSI4_PART1IMGF, i0,i1,i2)] + diagnostic_output_gfs[IDX4S(PSI4_PART2IMGF, i0,i1,i2)]); // 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; } } // Step 3: Perform integrations across all l,m modes from l=2 up to and including L_MAX (global variable): lowlevel_decompose_psi4_into_swm2_modes(Nxx_plus_2NGHOSTS1,Nxx_plus_2NGHOSTS2, dxx1,dxx2, time, R_ext, th_array, sinth_array, ph_array, psi4r_at_R_ext,psi4i_at_R_ext); } // Step 4: Free all allocated memory: free(psi4r_at_R_ext); free(psi4i_at_R_ext); free(sinth_array); free(th_array); free(ph_array); """ print_msg_with_timing("Spin-weight s=-2 Spherical Harmonics", msg="Ccodegen", startstop="stop", starttime=starttime) add_to_Cfunction_dict(includes=includes, prefunc=prefunc, desc=desc, name=name, params=params, body=body, rel_path_to_Cparams=rel_path_to_Cparams) return pickle_NRPy_env()
def add_psi4_part_to_Cfunction_dict(includes=None, rel_path_to_Cparams=os.path.join("."), whichpart=0, setPsi4tozero=False, OMP_pragma_on="i2"): starttime = print_msg_with_timing("psi4, part " + str(whichpart), msg="Ccodegen", startstop="start") # Set up the C function for psi4 if includes is None: includes = [] includes += ["NRPy_function_prototypes.h"] desc = "Compute psi4 at all interior gridpoints, part " + str(whichpart) name = "psi4_part" + str(whichpart) params = """const paramstruct *restrict params, const REAL *restrict in_gfs, REAL *restrict xx[3], REAL *restrict aux_gfs""" body = "" gri.register_gridfunctions("AUX", [ "psi4_part" + str(whichpart) + "re", "psi4_part" + str(whichpart) + "im" ]) FD_outCparams = "outCverbose=False,enable_SIMD=False,CSE_sorting=none" if not setPsi4tozero: # Set the body of the function # First compute the symbolic expressions psi4.Psi4(specify_tetrad=False) # We really don't want to store these "Cparameters" permanently; they'll be set via function call... # so we make a copy of the original par.glb_Cparams_list (sans tetrad vectors) and restore it below Cparams_list_orig = par.glb_Cparams_list.copy() par.Cparameters("REAL", __name__, ["mre4U0", "mre4U1", "mre4U2", "mre4U3"], [0, 0, 0, 0]) par.Cparameters("REAL", __name__, ["mim4U0", "mim4U1", "mim4U2", "mim4U3"], [0, 0, 0, 0]) par.Cparameters("REAL", __name__, ["n4U0", "n4U1", "n4U2", "n4U3"], [0, 0, 0, 0]) body += """ REAL mre4U0,mre4U1,mre4U2,mre4U3,mim4U0,mim4U1,mim4U2,mim4U3,n4U0,n4U1,n4U2,n4U3; psi4_tetrad(params, in_gfs[IDX4S(CFGF, i0,i1,i2)], in_gfs[IDX4S(HDD00GF, i0,i1,i2)], in_gfs[IDX4S(HDD01GF, i0,i1,i2)], in_gfs[IDX4S(HDD02GF, i0,i1,i2)], in_gfs[IDX4S(HDD11GF, i0,i1,i2)], in_gfs[IDX4S(HDD12GF, i0,i1,i2)], in_gfs[IDX4S(HDD22GF, i0,i1,i2)], &mre4U0,&mre4U1,&mre4U2,&mre4U3,&mim4U0,&mim4U1,&mim4U2,&mim4U3,&n4U0,&n4U1,&n4U2,&n4U3, xx, i0,i1,i2); """ body += "REAL xCart_rel_to_globalgrid_center[3];\n" body += "xx_to_Cart(params, xx, i0, i1, i2, xCart_rel_to_globalgrid_center);\n" body += "int ignore_Cart_to_i0i1i2[3]; REAL xx_rel_to_globalgridorigin[3];\n" body += "Cart_to_xx_and_nearest_i0i1i2_global_grid_center(params, xCart_rel_to_globalgrid_center,xx_rel_to_globalgridorigin,ignore_Cart_to_i0i1i2);\n" for i in range(3): body += "const REAL xx" + str( i) + "=xx_rel_to_globalgridorigin[" + str(i) + "];\n" body += fin.FD_outputC("returnstring", [ lhrh(lhs=gri.gfaccess("in_gfs", "psi4_part" + str(whichpart) + "re"), rhs=psi4.psi4_re_pt[whichpart]), lhrh(lhs=gri.gfaccess("in_gfs", "psi4_part" + str(whichpart) + "im"), rhs=psi4.psi4_im_pt[whichpart]) ], params=FD_outCparams) par.glb_Cparams_list = Cparams_list_orig.copy() elif setPsi4tozero: body += fin.FD_outputC("returnstring", [ lhrh(lhs=gri.gfaccess("in_gfs", "psi4_part" + str(whichpart) + "re"), rhs=sp.sympify(0)), lhrh(lhs=gri.gfaccess("in_gfs", "psi4_part" + str(whichpart) + "im"), rhs=sp.sympify(0)) ], params=FD_outCparams) enable_SIMD = False enable_rfm_precompute = False print_msg_with_timing("psi4, part " + str(whichpart), msg="Ccodegen", startstop="stop", starttime=starttime) add_to_Cfunction_dict(includes=includes, desc=desc, name=name, params=params, body=body, loopopts=get_loopopts("InteriorPoints", enable_SIMD, enable_rfm_precompute, OMP_pragma_on, enable_xxs=False), rel_path_to_Cparams=rel_path_to_Cparams) return pickle_NRPy_env()
def add_enforce_detgammahat_constraint_to_Cfunction_dict( includes=None, rel_path_to_Cparams=os.path.join("."), enable_rfm_precompute=True, enable_golden_kernels=False, OMP_pragma_on="i2", func_name_suffix=""): # This function disables SIMD, as it includes cbrt() and abs() functions. if includes is None: includes = [] # This function does not use finite differences! # enable_FD_functions = bool(par.parval_from_str("finite_difference::enable_FD_functions")) # if enable_FD_functions: # includes += ["finite_difference_functions.h"] # Set up the C function for enforcing the det(gammabar) = det(gammahat) BSSN algebraic constraint desc = "Enforce the det(gammabar) = det(gammahat) (algebraic) constraint" name = "enforce_detgammahat_constraint" + func_name_suffix params = "const paramstruct *restrict params, " if enable_rfm_precompute: params += "const rfm_struct *restrict rfmstruct, " else: params += "REAL *xx[3], " params += "REAL *restrict in_gfs" # Construct body: enforce_detg_constraint_symb_expressions = EGC.Enforce_Detgammahat_Constraint_symb_expressions( ) preloop = "" enableCparameters = True # Set up preloop in case we're outputting code for the Einstein Toolkit (ETK) if par.parval_from_str("grid::GridFuncMemAccess") == "ETK": params, preloop = set_ETK_func_params_preloop(func_name_suffix, enable_SIMD=False) enableCparameters = False FD_outCparams = "outCverbose=False,enable_SIMD=False" FD_outCparams += ",GoldenKernelsEnable=" + str(enable_golden_kernels) starttime = print_msg_with_timing( "Enforcing det(gammabar)=det(gammahat) constraint", msg="Ccodegen", startstop="start") body = fin.FD_outputC("returnstring", enforce_detg_constraint_symb_expressions, params=FD_outCparams) print_msg_with_timing("Enforcing det(gammabar)=det(gammahat) constraint", msg="Ccodegen", startstop="stop", starttime=starttime) enable_SIMD = False add_to_Cfunction_dict(includes=includes, desc=desc, name=name, params=params, preloop=preloop, body=body, loopopts=get_loopopts("AllPoints", enable_SIMD, enable_rfm_precompute, OMP_pragma_on), rel_path_to_Cparams=rel_path_to_Cparams, enableCparameters=enableCparameters) return pickle_NRPy_env()
def add_BSSN_constraints_to_Cfunction_dict( includes=None, rel_path_to_Cparams=os.path.join("."), enable_rfm_precompute=True, enable_golden_kernels=False, enable_SIMD=True, enable_stress_energy_source_terms=False, leave_Ricci_symbolic=True, output_H_only=False, OMP_pragma_on="i2", func_name_suffix=""): if includes is None: includes = [] if enable_SIMD: includes += [os.path.join("SIMD", "SIMD_intrinsics.h")] enable_FD_functions = bool( par.parval_from_str("finite_difference::enable_FD_functions")) if enable_FD_functions: includes += ["finite_difference_functions.h"] # Set up the C function for the BSSN constraints desc = "Evaluate the BSSN constraints" name = "BSSN_constraints" + func_name_suffix params = "const paramstruct *restrict params, " if enable_rfm_precompute: params += "const rfm_struct *restrict rfmstruct, " else: params += "REAL *xx[3], " params += """ const REAL *restrict in_gfs, const REAL *restrict auxevol_gfs, REAL *restrict aux_gfs""" # Construct body: BSSN_constraints_SymbExpressions = BSSN_constraints__generate_symbolic_expressions( enable_stress_energy_source_terms, leave_Ricci_symbolic=leave_Ricci_symbolic, output_H_only=output_H_only) preloop = "" enableCparameters = True # Set up preloop in case we're outputting code for the Einstein Toolkit (ETK) if par.parval_from_str("grid::GridFuncMemAccess") == "ETK": params, preloop = set_ETK_func_params_preloop(func_name_suffix) enableCparameters = False FD_outCparams = "outCverbose=False,enable_SIMD=" + str(enable_SIMD) FD_outCparams += ",GoldenKernelsEnable=" + str(enable_golden_kernels) FDorder = par.parval_from_str("finite_difference::FD_CENTDERIVS_ORDER") starttime = print_msg_with_timing("BSSN constraints (FD order=" + str(FDorder) + ")", msg="Ccodegen", startstop="start") body = fin.FD_outputC("returnstring", BSSN_constraints_SymbExpressions, params=FD_outCparams) print_msg_with_timing("BSSN constraints (FD order=" + str(FDorder) + ")", msg="Ccodegen", startstop="stop", starttime=starttime) add_to_Cfunction_dict(includes=includes, desc=desc, name=name, params=params, preloop=preloop, body=body, loopopts=get_loopopts("InteriorPoints", enable_SIMD, enable_rfm_precompute, OMP_pragma_on), rel_path_to_Cparams=rel_path_to_Cparams, enableCparameters=enableCparameters) return pickle_NRPy_env()
def add_Ricci_eval_to_Cfunction_dict( includes=None, rel_path_to_Cparams=os.path.join("."), enable_rfm_precompute=True, enable_golden_kernels=False, enable_SIMD=True, enable_split_for_optimizations_doesnt_help=False, OMP_pragma_on="i2", func_name_suffix=""): if includes is None: includes = [] if enable_SIMD: includes += [os.path.join("SIMD", "SIMD_intrinsics.h")] enable_FD_functions = bool( par.parval_from_str("finite_difference::enable_FD_functions")) if enable_FD_functions: includes += ["finite_difference_functions.h"] # Set up the C function for the 3-Ricci tensor desc = "Evaluate the 3-Ricci tensor" name = "Ricci_eval" + func_name_suffix params = "const paramstruct *restrict params, " if enable_rfm_precompute: params += "const rfm_struct *restrict rfmstruct, " else: params += "REAL *xx[3], " params += "const REAL *restrict in_gfs, REAL *restrict auxevol_gfs" # Construct body: Ricci_SymbExpressions = Ricci__generate_symbolic_expressions() FD_outCparams = "outCverbose=False,enable_SIMD=" + str(enable_SIMD) FD_outCparams += ",GoldenKernelsEnable=" + str(enable_golden_kernels) loopopts = get_loopopts("InteriorPoints", enable_SIMD, enable_rfm_precompute, OMP_pragma_on) FDorder = par.parval_from_str("finite_difference::FD_CENTDERIVS_ORDER") starttime = print_msg_with_timing("3-Ricci tensor (FD order=" + str(FDorder) + ")", msg="Ccodegen", startstop="start") # Construct body: preloop = "" enableCparameters = True # Set up preloop in case we're outputting code for the Einstein Toolkit (ETK) if par.parval_from_str("grid::GridFuncMemAccess") == "ETK": params, preloop = set_ETK_func_params_preloop(func_name_suffix) enableCparameters = False if enable_split_for_optimizations_doesnt_help and FDorder >= 8: loopopts += ",DisableOpenMP" Ricci_SymbExpressions_pt1 = [] Ricci_SymbExpressions_pt2 = [] for lhsrhs in Ricci_SymbExpressions: if "RBARDD00" in lhsrhs.lhs or "RBARDD11" in lhsrhs.lhs or "RBARDD22" in lhsrhs.lhs: Ricci_SymbExpressions_pt1.append( lhrh(lhs=lhsrhs.lhs, rhs=lhsrhs.rhs)) else: Ricci_SymbExpressions_pt2.append( lhrh(lhs=lhsrhs.lhs, rhs=lhsrhs.rhs)) preloop = """#pragma omp parallel { #pragma omp for """ preloopbody = fin.FD_outputC("returnstring", Ricci_SymbExpressions_pt1, params=FD_outCparams) preloop += lp.simple_loop(loopopts, preloopbody) preloop += "#pragma omp for\n" body = fin.FD_outputC("returnstring", Ricci_SymbExpressions_pt2, params=FD_outCparams) postloop = "\n } // END #pragma omp parallel\n" else: body = fin.FD_outputC("returnstring", Ricci_SymbExpressions, params=FD_outCparams) postloop = "" print_msg_with_timing("3-Ricci tensor (FD order=" + str(FDorder) + ")", msg="Ccodegen", startstop="stop", starttime=starttime) add_to_Cfunction_dict(includes=includes, desc=desc, name=name, params=params, preloop=preloop, body=body, loopopts=loopopts, postloop=postloop, rel_path_to_Cparams=rel_path_to_Cparams, enableCparameters=enableCparameters) return pickle_NRPy_env()
def add_rhs_eval_to_Cfunction_dict( includes=None, rel_path_to_Cparams=os.path.join("."), enable_rfm_precompute=True, enable_golden_kernels=False, enable_SIMD=True, enable_split_for_optimizations_doesnt_help=False, LapseCondition="OnePlusLog", ShiftCondition="GammaDriving2ndOrder_Covariant", enable_KreissOliger_dissipation=False, enable_stress_energy_source_terms=False, leave_Ricci_symbolic=True, OMP_pragma_on="i2", func_name_suffix=""): if includes is None: includes = [] if enable_SIMD: includes += [os.path.join("SIMD", "SIMD_intrinsics.h")] enable_FD_functions = bool( par.parval_from_str("finite_difference::enable_FD_functions")) if enable_FD_functions: includes += ["finite_difference_functions.h"] # Set up the C function for the BSSN RHSs desc = "Evaluate the BSSN RHSs" name = "rhs_eval" + func_name_suffix params = "const paramstruct *restrict params, " if enable_rfm_precompute: params += "const rfm_struct *restrict rfmstruct, " else: params += "REAL *xx[3], " params += """ const REAL *restrict auxevol_gfs,const REAL *restrict in_gfs,REAL *restrict rhs_gfs""" betaU, BSSN_RHSs_SymbExpressions = \ BSSN_RHSs__generate_symbolic_expressions(LapseCondition=LapseCondition, ShiftCondition=ShiftCondition, enable_KreissOliger_dissipation=enable_KreissOliger_dissipation, enable_stress_energy_source_terms=enable_stress_energy_source_terms, leave_Ricci_symbolic=leave_Ricci_symbolic) # Construct body: preloop = "" enableCparameters = True # Set up preloop in case we're outputting code for the Einstein Toolkit (ETK) if par.parval_from_str("grid::GridFuncMemAccess") == "ETK": params, preloop = set_ETK_func_params_preloop(func_name_suffix) enableCparameters = False FD_outCparams = "outCverbose=False,enable_SIMD=" + str(enable_SIMD) FD_outCparams += ",GoldenKernelsEnable=" + str(enable_golden_kernels) loopopts = get_loopopts("InteriorPoints", enable_SIMD, enable_rfm_precompute, OMP_pragma_on) FDorder = par.parval_from_str("finite_difference::FD_CENTDERIVS_ORDER") starttime = print_msg_with_timing("BSSN_RHSs (FD order=" + str(FDorder) + ")", msg="Ccodegen", startstop="start") if enable_split_for_optimizations_doesnt_help and FDorder == 6: loopopts += ",DisableOpenMP" BSSN_RHSs_SymbExpressions_pt1 = [] BSSN_RHSs_SymbExpressions_pt2 = [] for lhsrhs in BSSN_RHSs_SymbExpressions: if "BETU" in lhsrhs.lhs or "LAMBDAU" in lhsrhs.lhs: BSSN_RHSs_SymbExpressions_pt1.append( lhrh(lhs=lhsrhs.lhs, rhs=lhsrhs.rhs)) else: BSSN_RHSs_SymbExpressions_pt2.append( lhrh(lhs=lhsrhs.lhs, rhs=lhsrhs.rhs)) preloop += """#pragma omp parallel { """ preloopbody = fin.FD_outputC("returnstring", BSSN_RHSs_SymbExpressions_pt1, params=FD_outCparams, upwindcontrolvec=betaU) preloop += "\n#pragma omp for\n" + lp.simple_loop( loopopts, preloopbody) preloop += "\n#pragma omp for\n" body = fin.FD_outputC("returnstring", BSSN_RHSs_SymbExpressions_pt2, params=FD_outCparams, upwindcontrolvec=betaU) postloop = "\n } // END #pragma omp parallel\n" else: preloop += "" body = fin.FD_outputC("returnstring", BSSN_RHSs_SymbExpressions, params=FD_outCparams, upwindcontrolvec=betaU) postloop = "" print_msg_with_timing("BSSN_RHSs (FD order=" + str(FDorder) + ")", msg="Ccodegen", startstop="stop", starttime=starttime) add_to_Cfunction_dict(includes=includes, desc=desc, name=name, params=params, preloop=preloop, body=body, loopopts=loopopts, postloop=postloop, rel_path_to_Cparams=rel_path_to_Cparams, enableCparameters=enableCparameters) return pickle_NRPy_env()