Ejemplo n.º 1
0
    def Ricci__generate_symbolic_expressions():
        ######################################
        # START: GENERATE SYMBOLIC EXPRESSIONS
        print("Generating symbolic expressions for Ricci tensor...")
        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).
        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(outdir,"rfm_files/"))

        # Evaluate BSSN + BSSN gauge RHSs with rfm_precompute enabled:
        import BSSN.BSSN_quantities as Bq

        # Next compute Ricci tensor
        # par.set_parval_from_str("BSSN.BSSN_quantities::LeaveRicciSymbolic","False")
        RbarDD_already_registered = False
        for i in range(len(gri.glb_gridfcs_list)):
            if "RbarDD00" in gri.glb_gridfcs_list[i].name:
                RbarDD_already_registered = True
        if not RbarDD_already_registered:
            # We ignore the return value of ixp.register_gridfunctions_for_single_rank2() below
            #    as it is unused.
            ixp.register_gridfunctions_for_single_rank2("AUXEVOL","RbarDD","sym01")
        rhs.BSSN_RHSs()
        Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()

        # 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 Ricci symbolic expressions in "+str(end-start)+" seconds.")
        # END: GENERATE SYMBOLIC EXPRESSIONS
        ######################################

        Ricci_SymbExpressions = [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])]

        return [Ricci_SymbExpressions]
Ejemplo n.º 2
0
def Ricci__generate_symbolic_expressions():
    ######################################
    # START: GENERATE SYMBOLIC EXPRESSIONS
    starttime = print_msg_with_timing("3-Ricci tensor",
                                      msg="Symbolic",
                                      startstop="start")

    # Evaluate 3-Ricci tensor:
    import BSSN.BSSN_quantities as Bq
    par.set_parval_from_str("BSSN.BSSN_quantities::LeaveRicciSymbolic",
                            "False")

    # Register all BSSN gridfunctions if not registered already
    Bq.BSSN_basic_tensors()
    # Next compute Ricci tensor
    Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()
    # END: GENERATE SYMBOLIC EXPRESSIONS
    ######################################
    # Must register RbarDD as gridfunctions, as we're outputting them to gridfunctions here:
    foundit = False
    for i in range(len(gri.glb_gridfcs_list)):
        if "RbarDD00" in gri.glb_gridfcs_list[i].name:
            foundit = True
    if not foundit:
        ixp.register_gridfunctions_for_single_rank2("AUXEVOL", "RbarDD",
                                                    "sym01")

    Ricci_SymbExpressions = [
        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])
    ]
    print_msg_with_timing("3-Ricci tensor",
                          msg="Symbolic",
                          startstop="stop",
                          starttime=starttime)
    return Ricci_SymbExpressions
Ejemplo n.º 3
0
def BSSN_RHSs():
    # Step 1.c: Given the chosen coordinate system, set up
    #           corresponding reference metric and needed
    #           reference metric quantities
    # The following function call sets up the reference metric
    #    and related quantities, including rescaling matrices ReDD,
    #    ReU, and hatted quantities.
    rfm.reference_metric()

    global have_already_called_BSSN_RHSs_function  # setting to global enables other modules to see updated value.
    have_already_called_BSSN_RHSs_function = True

    # Step 1.d: Set spatial dimension (must be 3 for BSSN, as BSSN is
    #           a 3+1-dimensional decomposition of the general
    #           relativistic field equations)
    DIM = 3

    # Step 1.e: Import all basic (unrescaled) BSSN scalars & tensors
    import BSSN.BSSN_quantities as Bq
    Bq.BSSN_basic_tensors()
    gammabarDD = Bq.gammabarDD
    AbarDD = Bq.AbarDD
    LambdabarU = Bq.LambdabarU
    trK = Bq.trK
    alpha = Bq.alpha
    betaU = Bq.betaU

    # Step 1.f: Import all neeeded rescaled BSSN tensors:
    aDD = Bq.aDD
    cf = Bq.cf
    lambdaU = Bq.lambdaU

    # Step 2.a.i: Import derivative expressions for betaU defined in the BSSN.BSSN_quantities module:
    Bq.betaU_derivs()
    betaU_dD = Bq.betaU_dD
    betaU_dDD = Bq.betaU_dDD
    # Step 2.a.ii: Import derivative expression for gammabarDD
    Bq.gammabar__inverse_and_derivs()
    gammabarDD_dupD = Bq.gammabarDD_dupD

    # Step 2.a.iii: First term of \partial_t \bar{\gamma}_{i j} right-hand side:
    # \beta^k \bar{\gamma}_{ij,k} + \beta^k_{,i} \bar{\gamma}_{kj} + \beta^k_{,j} \bar{\gamma}_{ik}
    gammabar_rhsDD = ixp.zerorank2()
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                gammabar_rhsDD[i][j] += betaU[k] * gammabarDD_dupD[i][j][k] + betaU_dD[k][i] * gammabarDD[k][j] \
                                        + betaU_dD[k][j] * gammabarDD[i][k]

    # Step 2.b.i: First import \bar{A}_{ij} = AbarDD[i][j], and its contraction trAbar = \bar{A}^k_k
    #           from BSSN.BSSN_quantities
    Bq.AbarUU_AbarUD_trAbar_AbarDD_dD()
    trAbar = Bq.trAbar

    # Step 2.b.ii: Import detgammabar quantities from BSSN.BSSN_quantities:
    Bq.detgammabar_and_derivs()
    detgammabar = Bq.detgammabar
    detgammabar_dD = Bq.detgammabar_dD

    # Step 2.b.ii: Compute the contraction \bar{D}_k \beta^k = \beta^k_{,k} + \frac{\beta^k \bar{\gamma}_{,k}}{2 \bar{\gamma}}
    Dbarbetacontraction = sp.sympify(0)
    for k in range(DIM):
        Dbarbetacontraction += betaU_dD[k][
            k] + betaU[k] * detgammabar_dD[k] / (2 * detgammabar)

    # Step 2.b.iii: Second term of \partial_t \bar{\gamma}_{i j} right-hand side:
    # \frac{2}{3} \bar{\gamma}_{i j} \left (\alpha \bar{A}_{k}^{k} - \bar{D}_{k} \beta^{k}\right )
    for i in range(DIM):
        for j in range(DIM):
            gammabar_rhsDD[i][j] += sp.Rational(2, 3) * gammabarDD[i][j] * (
                alpha * trAbar - Dbarbetacontraction)

    # Step 2.c: Third term of \partial_t \bar{\gamma}_{i j} right-hand side:
    # -2 \alpha \bar{A}_{ij}
    for i in range(DIM):
        for j in range(DIM):
            gammabar_rhsDD[i][j] += -2 * alpha * AbarDD[i][j]

    # Step 3.a: First term of \partial_t \bar{A}_{i j}:
    # \beta^k \partial_k \bar{A}_{ij} + \partial_i \beta^k \bar{A}_{kj} + \partial_j \beta^k \bar{A}_{ik}

    # First define AbarDD_dupD:
    AbarDD_dupD = Bq.AbarDD_dupD  # From Bq.AbarUU_AbarUD_trAbar_AbarDD_dD()

    Abar_rhsDD = ixp.zerorank2()
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                Abar_rhsDD[i][j] += betaU[k] * AbarDD_dupD[i][j][k] + betaU_dD[k][i] * AbarDD[k][j] \
                                    + betaU_dD[k][j] * AbarDD[i][k]

    # Step 3.b: Second term of \partial_t \bar{A}_{i j}:
    # - (2/3) \bar{A}_{i j} \bar{D}_{k} \beta^{k} - 2 \alpha \bar{A}_{i k} {\bar{A}^{k}}_{j} + \alpha \bar{A}_{i j} K
    gammabarUU = Bq.gammabarUU  # From Bq.gammabar__inverse_and_derivs()
    AbarUD = Bq.AbarUD  # From Bq.AbarUU_AbarUD_trAbar()
    for i in range(DIM):
        for j in range(DIM):
            Abar_rhsDD[i][j] += -sp.Rational(2, 3) * AbarDD[i][
                j] * Dbarbetacontraction + alpha * AbarDD[i][j] * trK
            for k in range(DIM):
                Abar_rhsDD[i][j] += -2 * alpha * AbarDD[i][k] * AbarUD[k][j]

    # Step 3.c.i: Define partial derivatives of \phi in terms of evolved quantity "cf":
    Bq.phi_and_derivs()
    phi_dD = Bq.phi_dD
    phi_dupD = Bq.phi_dupD
    phi_dDD = Bq.phi_dDD
    exp_m4phi = Bq.exp_m4phi
    phi_dBarD = Bq.phi_dBarD  # phi_dBarD = Dbar_i phi = phi_dD (since phi is a scalar)
    phi_dBarDD = Bq.phi_dBarDD  # phi_dBarDD = Dbar_i Dbar_j phi (covariant derivative)

    # Step 3.c.ii: Define RbarDD
    Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()
    RbarDD = Bq.RbarDD

    # Step 3.c.iii: Define first and second derivatives of \alpha, as well as
    #         \bar{D}_i \bar{D}_j \alpha, which is defined just like phi
    alpha_dD = ixp.declarerank1("alpha_dD")
    alpha_dDD = ixp.declarerank2("alpha_dDD", "sym01")
    alpha_dBarD = alpha_dD
    alpha_dBarDD = ixp.zerorank2()
    GammabarUDD = Bq.GammabarUDD  # Defined in Bq.gammabar__inverse_and_derivs()
    for i in range(DIM):
        for j in range(DIM):
            alpha_dBarDD[i][j] = alpha_dDD[i][j]
            for k in range(DIM):
                alpha_dBarDD[i][j] += -GammabarUDD[k][i][j] * alpha_dD[k]

    # Step 3.c.iv: Define the terms in curly braces:
    curlybrackettermsDD = ixp.zerorank2()
    for i in range(DIM):
        for j in range(DIM):
            curlybrackettermsDD[i][j] = -2 * alpha * phi_dBarDD[i][j] + 4 * alpha * phi_dBarD[i] * phi_dBarD[j] \
                                        + 2 * alpha_dBarD[i] * phi_dBarD[j] \
                                        + 2 * alpha_dBarD[j] * phi_dBarD[i] \
                                        - alpha_dBarDD[i][j] + alpha * RbarDD[i][j]

    # Step 3.c.v: Compute the trace:
    curlybracketterms_trace = sp.sympify(0)
    for i in range(DIM):
        for j in range(DIM):
            curlybracketterms_trace += gammabarUU[i][j] * curlybrackettermsDD[
                i][j]

    # Step 3.c.vi: Third and final term of Abar_rhsDD[i][j]:
    for i in range(DIM):
        for j in range(DIM):
            Abar_rhsDD[i][j] += exp_m4phi * (
                curlybrackettermsDD[i][j] -
                sp.Rational(1, 3) * gammabarDD[i][j] * curlybracketterms_trace)

    # Step 4: Right-hand side of conformal factor variable "cf". Supported
    #          options include: cf=phi, cf=W=e^(-2*phi) (default), and cf=chi=e^(-4*phi)
    # \partial_t phi = \left[\beta^k \partial_k \phi \right] <- TERM 1
    #                  + \frac{1}{6} \left (\bar{D}_{k} \beta^{k} - \alpha K \right ) <- TERM 2
    global cf_rhs
    cf_rhs = sp.Rational(1, 6) * (Dbarbetacontraction - alpha * trK)  # Term 2
    for k in range(DIM):
        cf_rhs += betaU[k] * phi_dupD[k]  # Term 1

    # Next multiply to convert phi_rhs to cf_rhs.
    if par.parval_from_str(
            "BSSN.BSSN_quantities::EvolvedConformalFactor_cf") == "phi":
        pass  # do nothing; cf_rhs = phi_rhs
    elif par.parval_from_str(
            "BSSN.BSSN_quantities::EvolvedConformalFactor_cf") == "W":
        cf_rhs *= -2 * cf  # cf_rhs = -2*cf*phi_rhs
    elif par.parval_from_str(
            "BSSN.BSSN_quantities::EvolvedConformalFactor_cf") == "chi":
        cf_rhs *= -4 * cf  # cf_rhs = -4*cf*phi_rhs
    else:
        print("Error: EvolvedConformalFactor_cf == " + par.parval_from_str(
            "BSSN.BSSN_quantities::EvolvedConformalFactor_cf") +
              " unsupported!")
        exit(1)

    # Step 5: right-hand side of trK (trace of extrinsic curvature):
    # \partial_t K = \beta^k \partial_k K <- TERM 1
    #           + \frac{1}{3} \alpha K^{2} <- TERM 2
    #           + \alpha \bar{A}_{i j} \bar{A}^{i j} <- TERM 3
    #           - - e^{-4 \phi} (\bar{D}_{i} \bar{D}^{i} \alpha + 2 \bar{D}^{i} \alpha \bar{D}_{i} \phi ) <- TERM 4
    global trK_rhs
    # TERM 2:
    trK_rhs = sp.Rational(1, 3) * alpha * trK * trK
    trK_dupD = ixp.declarerank1("trK_dupD")
    for i in range(DIM):
        # TERM 1:
        trK_rhs += betaU[i] * trK_dupD[i]
    for i in range(DIM):
        for j in range(DIM):
            # TERM 4:
            trK_rhs += -exp_m4phi * gammabarUU[i][j] * (
                alpha_dBarDD[i][j] + 2 * alpha_dBarD[j] * phi_dBarD[i])
    AbarUU = Bq.AbarUU  # From Bq.AbarUU_AbarUD_trAbar()
    for i in range(DIM):
        for j in range(DIM):
            # TERM 3:
            trK_rhs += alpha * AbarDD[i][j] * AbarUU[i][j]

    # Step 6: right-hand side of \partial_t \bar{\Lambda}^i:
    # \partial_t \bar{\Lambda}^i = \beta^k \partial_k \bar{\Lambda}^i - \partial_k \beta^i \bar{\Lambda}^k <- TERM 1
    #                            + \bar{\gamma}^{j k} \hat{D}_{j} \hat{D}_{k} \beta^{i} <- TERM 2
    #                            + \frac{2}{3} \Delta^{i} \bar{D}_{j} \beta^{j} <- TERM 3
    #                            + \frac{1}{3} \bar{D}^{i} \bar{D}_{j} \beta^{j} <- TERM 4
    #                            - 2 \bar{A}^{i j} (\partial_{j} \alpha - 6 \partial_{j} \phi) <- TERM 5
    #                            + 2 \alpha \bar{A}^{j k} \Delta_{j k}^{i} <- TERM 6
    #                            - \frac{4}{3} \alpha \bar{\gamma}^{i j} \partial_{j} K <- TERM 7

    # Step 6.a: Term 1 of \partial_t \bar{\Lambda}^i: \beta^k \partial_k \bar{\Lambda}^i - \partial_k \beta^i \bar{\Lambda}^k
    # First we declare \bar{\Lambda}^i and \bar{\Lambda}^i_{,j} in terms of \lambda^i and \lambda^i_{,j}
    global LambdabarU_dupD  # Used on the RHS of the Gamma-driving shift conditions
    LambdabarU_dupD = ixp.zerorank2()
    lambdaU_dupD = ixp.declarerank2("lambdaU_dupD", "nosym")
    for i in range(DIM):
        for j in range(DIM):
            LambdabarU_dupD[i][j] = lambdaU_dupD[i][j] * rfm.ReU[i] + lambdaU[
                i] * rfm.ReUdD[i][j]

    global Lambdabar_rhsU  # Used on the RHS of the Gamma-driving shift conditions
    Lambdabar_rhsU = ixp.zerorank1()
    for i in range(DIM):
        for k in range(DIM):
            Lambdabar_rhsU[i] += betaU[k] * LambdabarU_dupD[i][k] - betaU_dD[
                i][k] * LambdabarU[k]  # Term 1

    # Step 6.b: Term 2 of \partial_t \bar{\Lambda}^i = \bar{\gamma}^{jk} (Term 2a + Term 2b + Term 2c)
    # Term 2a: \bar{\gamma}^{jk} \beta^i_{,kj}
    Term2aUDD = ixp.zerorank3()
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                Term2aUDD[i][j][k] += betaU_dDD[i][k][j]
    # Term 2b: \hat{\Gamma}^i_{mk,j} \beta^m + \hat{\Gamma}^i_{mk} \beta^m_{,j}
    #          + \hat{\Gamma}^i_{dj}\beta^d_{,k} - \hat{\Gamma}^d_{kj} \beta^i_{,d}
    Term2bUDD = ixp.zerorank3()
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                for m in range(DIM):
                    Term2bUDD[i][j][k] += rfm.GammahatUDDdD[i][m][k][j] * betaU[m] \
                                          + rfm.GammahatUDD[i][m][k] * betaU_dD[m][j] \
                                          + rfm.GammahatUDD[i][m][j] * betaU_dD[m][k] \
                                          - rfm.GammahatUDD[m][k][j] * betaU_dD[i][m]
    # Term 2c: \hat{\Gamma}^i_{dj}\hat{\Gamma}^d_{mk} \beta^m - \hat{\Gamma}^d_{kj} \hat{\Gamma}^i_{md} \beta^m
    Term2cUDD = ixp.zerorank3()
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                for m in range(DIM):
                    for d in range(DIM):
                        Term2cUDD[i][j][k] += (rfm.GammahatUDD[i][d][j] * rfm.GammahatUDD[d][m][k] \
                                               - rfm.GammahatUDD[d][k][j] * rfm.GammahatUDD[i][m][d]) * betaU[m]

    Lambdabar_rhsUpieceU = ixp.zerorank1()

    # Put it all together to get Term 2:
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                Lambdabar_rhsU[i] += gammabarUU[j][k] * (Term2aUDD[i][j][k] +
                                                         Term2bUDD[i][j][k] +
                                                         Term2cUDD[i][j][k])
                Lambdabar_rhsUpieceU[i] += gammabarUU[j][k] * (
                    Term2aUDD[i][j][k] + Term2bUDD[i][j][k] +
                    Term2cUDD[i][j][k])

    # Step 6.c: Term 3 of \partial_t \bar{\Lambda}^i:
    #    \frac{2}{3} \Delta^{i} \bar{D}_{j} \beta^{j}
    DGammaU = Bq.DGammaU  # From Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()
    for i in range(DIM):
        Lambdabar_rhsU[i] += sp.Rational(
            2, 3) * DGammaU[i] * Dbarbetacontraction  # Term 3

    # Step 6.d: Term 4 of \partial_t \bar{\Lambda}^i:
    #           \frac{1}{3} \bar{D}^{i} \bar{D}_{j} \beta^{j}
    detgammabar_dDD = Bq.detgammabar_dDD  # From Bq.detgammabar_and_derivs()
    Dbarbetacontraction_dBarD = ixp.zerorank1()
    for k in range(DIM):
        for m in range(DIM):
            Dbarbetacontraction_dBarD[m] += betaU_dDD[k][k][m] + \
                                            (betaU_dD[k][m] * detgammabar_dD[k] +
                                             betaU[k] * detgammabar_dDD[k][m]) / (2 * detgammabar) \
                                            - betaU[k] * detgammabar_dD[k] * detgammabar_dD[m] / (
                                                        2 * detgammabar * detgammabar)
    for i in range(DIM):
        for m in range(DIM):
            Lambdabar_rhsU[i] += sp.Rational(
                1, 3) * gammabarUU[i][m] * Dbarbetacontraction_dBarD[m]

    # Step 6.e: Term 5 of \partial_t \bar{\Lambda}^i:
    #           - 2 \bar{A}^{i j} (\partial_{j} \alpha - 6 \alpha \partial_{j} \phi)
    for i in range(DIM):
        for j in range(DIM):
            Lambdabar_rhsU[i] += -2 * AbarUU[i][j] * (alpha_dD[j] -
                                                      6 * alpha * phi_dD[j])

    # Step 6.f: Term 6 of \partial_t \bar{\Lambda}^i:
    #           2 \alpha \bar{A}^{j k} \Delta^{i}_{j k}
    DGammaUDD = Bq.DGammaUDD  # From RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                Lambdabar_rhsU[
                    i] += 2 * alpha * AbarUU[j][k] * DGammaUDD[i][j][k]

    # Step 6.g: Term 7 of \partial_t \bar{\Lambda}^i:
    #           -\frac{4}{3} \alpha \bar{\gamma}^{i j} \partial_{j} K
    trK_dD = ixp.declarerank1("trK_dD")
    for i in range(DIM):
        for j in range(DIM):
            Lambdabar_rhsU[i] += -sp.Rational(
                4, 3) * alpha * gammabarUU[i][j] * trK_dD[j]

    # Step 7: Rescale the RHS quantities so that the evolved
    #         variables are smooth across coord singularities
    global h_rhsDD, a_rhsDD, lambda_rhsU
    h_rhsDD = ixp.zerorank2()
    a_rhsDD = ixp.zerorank2()
    lambda_rhsU = ixp.zerorank1()
    for i in range(DIM):
        lambda_rhsU[i] = Lambdabar_rhsU[i] / rfm.ReU[i]
        for j in range(DIM):
            h_rhsDD[i][j] = gammabar_rhsDD[i][j] / rfm.ReDD[i][j]
            a_rhsDD[i][j] = Abar_rhsDD[i][j] / rfm.ReDD[i][j]
Ejemplo n.º 4
0
    def BSSN_RHSs_Ricci__generate_symbolic_expressions():
        ######################################
        # START: GENERATE SYMBOLIC EXPRESSIONS
        # Store original finite-differencing order:
        FD_order_orig = par.parval_from_str("finite_difference::FD_CENTDERIVS_ORDER")
        # Set new finite-differencing order:
        par.set_parval_from_str("finite_difference::FD_CENTDERIVS_ORDER", FD_order)


        print("Generating symbolic expressions for BSSN RHSs and Ricci tensor...")
        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(outdir,"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(outdir,"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()

        if T4UU != None:
            import BSSN.BSSN_stress_energy_source_terms as Bsest
            Bsest.BSSN_source_terms_for_BSSN_RHSs(T4UU)
            rhs.trK_rhs += Bsest.sourceterm_trK_rhs
            for i in range(3):
                # Needed for Gamma-driving shift RHSs:
                rhs.Lambdabar_rhsU[i] += Bsest.sourceterm_Lambdabar_rhsU[i]
                # Needed for BSSN RHSs:
                rhs.lambda_rhsU[i]    += Bsest.sourceterm_lambda_rhsU[i]
                for j in range(3):
                    rhs.a_rhsDD[i][j] += Bsest.sourceterm_a_rhsDD[i][j]

        gaugerhs.BSSN_gauge_RHSs()

        # Add Kreiss-Oliger dissipation to the BSSN RHSs:
        thismodule = "KO_Dissipation"
        diss_strength = par.Cparameters("REAL", thismodule, "diss_strength", default_KO_strength)

        alpha_dKOD = ixp.declarerank1("alpha_dKOD")
        cf_dKOD    = ixp.declarerank1("cf_dKOD")
        trK_dKOD   = ixp.declarerank1("trK_dKOD")
        betU_dKOD    = ixp.declarerank2("betU_dKOD","nosym")
        vetU_dKOD    = ixp.declarerank2("vetU_dKOD","nosym")
        lambdaU_dKOD = ixp.declarerank2("lambdaU_dKOD","nosym")
        aDD_dKOD = ixp.declarerank3("aDD_dKOD","sym01")
        hDD_dKOD = ixp.declarerank3("hDD_dKOD","sym01")
        for k in range(3):
            gaugerhs.alpha_rhs += diss_strength*alpha_dKOD[k]*rfm.ReU[k] # ReU[k] = 1/scalefactor_orthog_funcform[k]
            rhs.cf_rhs         += diss_strength*   cf_dKOD[k]*rfm.ReU[k] # ReU[k] = 1/scalefactor_orthog_funcform[k]
            rhs.trK_rhs        += diss_strength*  trK_dKOD[k]*rfm.ReU[k] # ReU[k] = 1/scalefactor_orthog_funcform[k]
            for i in range(3):
                if "2ndOrder" in ShiftCondition:
                    gaugerhs.bet_rhsU[i] += diss_strength*   betU_dKOD[i][k]*rfm.ReU[k] # ReU[k] = 1/scalefactor_orthog_funcform[k]
                gaugerhs.vet_rhsU[i]     += diss_strength*   vetU_dKOD[i][k]*rfm.ReU[k] # ReU[k] = 1/scalefactor_orthog_funcform[k]
                rhs.lambda_rhsU[i]       += diss_strength*lambdaU_dKOD[i][k]*rfm.ReU[k] # ReU[k] = 1/scalefactor_orthog_funcform[k]
                for j in range(3):
                    rhs.a_rhsDD[i][j] += diss_strength*aDD_dKOD[i][j][k]*rfm.ReU[k] # ReU[k] = 1/scalefactor_orthog_funcform[k]
                    rhs.h_rhsDD[i][j] += diss_strength*hDD_dKOD[i][j][k]*rfm.ReU[k] # ReU[k] = 1/scalefactor_orthog_funcform[k]

        # We use betaU as our upwinding control vector:
        Bq.BSSN_basic_tensors()
        betaU = Bq.betaU

        # Next compute Ricci tensor
        par.set_parval_from_str("BSSN.BSSN_quantities::LeaveRicciSymbolic","False")
        Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()

        # 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("Finished BSSN symbolic expressions in "+str(end-start)+" seconds.")
        # Restore original finite-differencing order:
        par.set_parval_from_str("finite_difference::FD_CENTDERIVS_ORDER", FD_order)
        # END: GENERATE SYMBOLIC EXPRESSIONS
        ######################################

        BSSN_RHSs_SymbExpressions = [lhrh(lhs=gri.gfaccess("rhs_gfs","aDD00"),   rhs=rhs.a_rhsDD[0][0]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","aDD01"),   rhs=rhs.a_rhsDD[0][1]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","aDD02"),   rhs=rhs.a_rhsDD[0][2]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","aDD11"),   rhs=rhs.a_rhsDD[1][1]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","aDD12"),   rhs=rhs.a_rhsDD[1][2]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","aDD22"),   rhs=rhs.a_rhsDD[2][2]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","alpha"),   rhs=gaugerhs.alpha_rhs),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","betU0"),   rhs=gaugerhs.bet_rhsU[0]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","betU1"),   rhs=gaugerhs.bet_rhsU[1]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","betU2"),   rhs=gaugerhs.bet_rhsU[2]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","cf"),      rhs=rhs.cf_rhs),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","hDD00"),   rhs=rhs.h_rhsDD[0][0]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","hDD01")   ,rhs=rhs.h_rhsDD[0][1]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","hDD02"),   rhs=rhs.h_rhsDD[0][2]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","hDD11"),   rhs=rhs.h_rhsDD[1][1]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","hDD12"),   rhs=rhs.h_rhsDD[1][2]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","hDD22"),   rhs=rhs.h_rhsDD[2][2]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","lambdaU0"),rhs=rhs.lambda_rhsU[0]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","lambdaU1"),rhs=rhs.lambda_rhsU[1]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","lambdaU2"),rhs=rhs.lambda_rhsU[2]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","trK"),     rhs=rhs.trK_rhs),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","vetU0"),   rhs=gaugerhs.vet_rhsU[0]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","vetU1"),   rhs=gaugerhs.vet_rhsU[1]),
                                     lhrh(lhs=gri.gfaccess("rhs_gfs","vetU2"),   rhs=gaugerhs.vet_rhsU[2]) ]

        Ricci_SymbExpressions = [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])]

        return [betaU,BSSN_RHSs_SymbExpressions,Ricci_SymbExpressions]
Ejemplo n.º 5
0
def BSSN_constraints(add_T4UUmunu_source_terms=False):
    # Step 1.a: Set spatial dimension (must be 3 for BSSN, as BSSN is
    #           a 3+1-dimensional decomposition of the general
    #           relativistic field equations)
    DIM = 3

    # Step 1.b: Given the chosen coordinate system, set up
    #           corresponding reference metric and needed
    #           reference metric quantities
    # The following function call sets up the reference metric
    #    and related quantities, including rescaling matrices ReDD,
    #    ReU, and hatted quantities.
    rfm.reference_metric()

    # Step 2: Hamiltonian constraint.
    # First declare all needed variables
    Bq.declare_BSSN_gridfunctions_if_not_declared_already()  # Sets trK
    Bq.BSSN_basic_tensors()  # Sets AbarDD
    Bq.gammabar__inverse_and_derivs()  # Sets gammabarUU
    Bq.AbarUU_AbarUD_trAbar_AbarDD_dD()  # Sets AbarUU and AbarDD_dD
    Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()  # Sets RbarDD
    Bq.phi_and_derivs()  # Sets phi_dBarD & phi_dBarDD

    ###############################
    ###############################
    #  HAMILTONIAN CONSTRAINT
    ###############################
    ###############################

    # Term 1: 2/3 K^2
    global H
    H = sp.Rational(2, 3) * Bq.trK**2

    # Term 2: -A_{ij} A^{ij}
    for i in range(DIM):
        for j in range(DIM):
            H += -Bq.AbarDD[i][j] * Bq.AbarUU[i][j]

    # Term 3a: trace(Rbar)
    Rbartrace = sp.sympify(0)
    for i in range(DIM):
        for j in range(DIM):
            Rbartrace += Bq.gammabarUU[i][j] * Bq.RbarDD[i][j]

    # Term 3b: -8 \bar{\gamma}^{ij} \bar{D}_i \phi \bar{D}_j \phi = -8*phi_dBar_times_phi_dBar
    # Term 3c: -8 \bar{\gamma}^{ij} \bar{D}_i \bar{D}_j \phi      = -8*phi_dBarDD_contraction
    phi_dBar_times_phi_dBar = sp.sympify(0)  # Term 3b
    phi_dBarDD_contraction = sp.sympify(0)  # Term 3c
    for i in range(DIM):
        for j in range(DIM):
            phi_dBar_times_phi_dBar += Bq.gammabarUU[i][j] * Bq.phi_dBarD[
                i] * Bq.phi_dBarD[j]
            phi_dBarDD_contraction += Bq.gammabarUU[i][j] * Bq.phi_dBarDD[i][j]

    # Add Term 3:
    H += Bq.exp_m4phi * (Rbartrace - 8 *
                         (phi_dBar_times_phi_dBar + phi_dBarDD_contraction))

    if add_T4UUmunu_source_terms:
        M_PI = par.Cparameters("#define", thismodule, "M_PI",
                               "")  # M_PI is pi as defined in C
        BTmunu.define_BSSN_T4UUmunu_rescaled_source_terms()
        rho = BTmunu.rho
        H += -16 * M_PI * rho

    # FIXME: ADD T4UUmunu SOURCE TERMS TO MOMENTUM CONSTRAINT!

    # Step 3: M^i, the momentum constraint

    ###############################
    ###############################
    #  MOMENTUM CONSTRAINT
    ###############################
    ###############################

    # SEE Tutorial-BSSN_constraints.ipynb for full documentation.
    global MU
    MU = ixp.zerorank1()

    # Term 2: 6 A^{ij} \partial_j \phi:
    for i in range(DIM):
        for j in range(DIM):
            MU[i] += 6 * Bq.AbarUU[i][j] * Bq.phi_dD[j]

    # Term 3: -2/3 \bar{\gamma}^{ij} K_{,j}
    trK_dD = ixp.declarerank1(
        "trK_dD")  # Not defined in BSSN_RHSs; only trK_dupD is defined there.
    for i in range(DIM):
        for j in range(DIM):
            MU[i] += -sp.Rational(2, 3) * Bq.gammabarUU[i][j] * trK_dD[j]

    # First define aDD_dD:
    aDD_dD = ixp.declarerank3("aDD_dD", "sym01")

    # Then evaluate the conformal covariant derivative \bar{D}_j \bar{A}_{lm}
    AbarDD_dBarD = ixp.zerorank3()
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                AbarDD_dBarD[i][j][k] = Bq.AbarDD_dD[i][j][k]
                for l in range(DIM):
                    AbarDD_dBarD[i][j][
                        k] += -Bq.GammabarUDD[l][k][i] * Bq.AbarDD[l][j]
                    AbarDD_dBarD[i][j][
                        k] += -Bq.GammabarUDD[l][k][j] * Bq.AbarDD[i][l]

    # Term 1: Contract twice with the metric to make \bar{D}_{j} \bar{A}^{ij}
    for i in range(DIM):
        for j in range(DIM):
            for k in range(DIM):
                for l in range(DIM):
                    MU[i] += Bq.gammabarUU[i][k] * Bq.gammabarUU[j][
                        l] * AbarDD_dBarD[k][l][j]

    # Finally, we multiply by e^{-4 phi} and rescale the momentum constraint:
    for i in range(DIM):
        MU[i] *= Bq.exp_m4phi / rfm.ReU[i]
Ejemplo n.º 6
0
    def test_example_BSSN():
        parse_latex(r"""
            \begin{align}
                % keydef basis [x, y, z]
                % ignore "\\%", "\qquad"

                % vardef -kron 'deltaDD'
                % parse \hat{\gamma}_{ij} = \delta_{ij}
                % assign -diff_type=symbolic -metric 'gammahatDD'
                % vardef -diff_type=dD -symmetry=sym01 'hDD'
                % parse \bar{\gamma}_{ij} = h_{ij} + \hat{\gamma}_{ij}
                % assign -diff_type=dD -metric 'gammabarDD'

                % srepl "\beta" -> "\text{vet}"
                % vardef -diff_type=dD 'vetU'
                %% upwind pattern inside Lie derivative expansion
                % srepl -persist "\text{vet}^{<1>} \partial_{<1>}" -> "\text{vet}^{<1>} \vphantom{dupD} \partial_{<1>}"
                %% substitute tensor identity (see appropriate BSSN notebook)
                % srepl "\bar{D}_k \text{vet}^k" -> "(\partial_k \text{vet}^k + \frac{\partial_k \text{gammahatdet} \text{vet}^k}{2 \text{gammahatdet}})"

                % srepl "\bar{A}" -> "\text{a}"
                % vardef -diff_type=dD -symmetry=sym01 'aDD'
                % assign -metric='gammabar' 'aDD'
                % srepl "\partial_t \bar{\gamma}" -> "\text{h_rhs}"
                \partial_t \bar{\gamma}_{ij} &= \mathcal{L}_\beta \bar{\gamma}_{ij}
                    + \frac{2}{3} \bar{\gamma}_{ij} \left(\alpha \bar{A}^k{}_k - \bar{D}_k \beta^k\right) - 2 \alpha \bar{A}_{ij} \\

                % srepl "K" -> "\text{trK}"
                % vardef -diff_type=dD 'cf', 'trK'
                %% replace 'phi' with conformal factor cf = W = e^{-2\phi}
                % srepl "e^{-4\phi}" -> "\text{cf}^2"
                % srepl "\partial_t \phi = <1..> \\" -> "\text{cf_rhs} = -2 \text{cf} (<1..>) \\"
                % srepl -persist "\partial_{<1>} \phi" -> "\partial_{<1>} \text{cf} \frac{-1}{2 \text{cf}}"
                % srepl "\partial_<1> \phi" -> "\partial_<1> \text{cf} \frac{-1}{2 \text{cf}}"
                \partial_t \phi &= \mathcal{L}_\beta \phi + \frac{1}{6} \left(\bar{D}_k \beta^k - \alpha K \right) \\

                % vardef -diff_type=dD 'alpha'
                % srepl "\partial_t \text{trK}" -> "\text{trK_rhs}"
                \partial_t K &= \mathcal{L}_\beta K + \frac{1}{3} \alpha K^2 + \alpha \bar{A}_{ij} \bar{A}^{ij}
                    - e^{-4\phi} \left(\bar{D}_i \bar{D}^i \alpha + 2 \bar{D}^i \alpha \bar{D}_i \phi\right) \\

                % srepl "\bar{\Lambda}" -> "\text{lambda}"
                % vardef -diff_type=dD 'lambdaU'
                % parse \Delta^k_{ij} = \bar{\Gamma}^k_{ij} - \hat{\Gamma}^k_{ij}
                % assign -metric='gammabar' 'DeltaUDD'
                % parse \Delta^k = \bar{\gamma}^{ij} \Delta^k_{ij}
                % srepl "\partial_t \text{lambda}" -> "\text{Lambdabar_rhs}"
                \partial_t \bar{\Lambda}^i &= \mathcal{L}_\beta \bar{\Lambda}^i + \bar{\gamma}^{jk} \hat{D}_j \hat{D}_k \beta^i
                    + \frac{2}{3} \Delta^i \bar{D}_k \beta^k + \frac{1}{3} \bar{D}^i \bar{D}_k \beta^k \\%
                    &\qquad- 2 \bar{A}^{ij} \left(\partial_j \alpha - 6 \alpha \partial_j \phi\right)
                    + 2 \alpha \bar{A}^{jk} \Delta^i_{jk} - \frac{4}{3} \alpha \bar{\gamma}^{ij} \partial_j K \\

                % vardef -diff_type=dD -symmetry=sym01 'RbarDD'
                X_{ij} &= -2 \alpha \bar{D}_i \bar{D}_j \phi + 4 \alpha \bar{D}_i \phi \bar{D}_j \phi
                    + 2 \bar{D}_i \alpha \bar{D}_j \phi + 2 \bar{D}_j \alpha \bar{D}_i \phi
                    - \bar{D}_i \bar{D}_j \alpha + \alpha \bar{R}_{ij} \\
                \hat{X}_{ij} &= X_{ij} - \frac{1}{3} \bar{\gamma}_{ij} \bar{\gamma}^{kl} X_{kl} \\
                % srepl "\partial_t \text{a}" -> "\text{a_rhs}"
                \partial_t \bar{A}_{ij} &= \mathcal{L}_\beta \bar{A}_{ij} - \frac{2}{3} \bar{A}_{ij} \bar{D}_k \beta^k
                    - 2 \alpha \bar{A}_{ik} \bar{A}^k_j + \alpha \bar{A}_{ij} K + e^{-4\phi} \hat{X}_{ij} \\

                % srepl "\partial_t \alpha" -> "\text{alpha_rhs}"
                \partial_t \alpha &= \mathcal{L}_\beta \alpha - 2 \alpha K \\

                % srepl "B" -> "\text{bet}"
                % vardef -diff_type=dD 'betU'
                % srepl "\partial_t \text{vet}" -> "\text{vet_rhs}"
                \partial_t \beta^i &= \left[\beta^j \vphantom{dupD} \bar{D}_j \beta^i\right] + B^i \\

                % vardef -const 'eta'
                % srepl "\partial_t \text{bet}" -> "\text{bet_rhs}"
                \partial_t B^i &= \left[\beta^j \vphantom{dupD} \bar{D}_j B^i\right]
                    + \frac{3}{4} \left(\partial_t \bar{\Lambda}^i - \left[\beta^j \vphantom{dupD} \bar{D}_j \bar{\Lambda}^i\right]\right) - \eta B^i \\

                % parse \bar{R} = \bar{\gamma}^{ij} \bar{R}_{ij}
                % srepl "\bar{D}^2" -> "\bar{D}^i \bar{D}_i", "\mathcal{<1>}" -> "<1>"
                \mathcal{H} &= \frac{2}{3} K^2 - \bar{A}_{ij} \bar{A}^{ij}
                    + e^{-4\phi} \left(\bar{R} - 8 \bar{D}^i \phi \bar{D}_i \phi - 8 \bar{D}^2 \phi\right) \\

                \mathcal{M}^i &= e^{-4\phi} \left(\bar{D}_j \bar{A}^{ij} + 6 \bar{A}^{ij} \partial_j \phi
                    - \frac{2}{3} \bar{\gamma}^{ij} \partial_j K\right) \\

                \bar{R}_{ij} &= -\frac{1}{2} \bar{\gamma}^{kl} \hat{D}_k \hat{D}_l \bar{\gamma}_{ij}
                    + \frac{1}{2} \left(\bar{\gamma}_{ki} \hat{D}_j \bar{\Lambda}^k + \bar{\gamma}_{kj} \hat{D}_i \bar{\Lambda}^k\right)
                    + \frac{1}{2} \Delta^k \left(\Delta_{ijk} + \Delta_{jik}\right) \\%
                    &\qquad+ \bar{\gamma}^{kl} \left(\Delta^m_{ki} \Delta_{jml} + \Delta^m_{kj} \Delta_{iml} + \Delta^m_{ik} \Delta_{mjl}\right)
            \end{align}
        """,
                    ignore_warning=True)
        par.set_parval_from_str('reference_metric::CoordSystem', 'Cartesian')
        par.set_parval_from_str('BSSN.BSSN_quantities::LeaveRicciSymbolic',
                                'True')
        rfm.reference_metric()
        Brhs.BSSN_RHSs()
        gaugerhs.BSSN_gauge_RHSs()
        bssncon.BSSN_constraints()
        par.set_parval_from_str('BSSN.BSSN_quantities::LeaveRicciSymbolic',
                                'False')
        Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()
        assert_equal(
            {
                'h_rhsDD': h_rhsDD,
                'cf_rhs': cf_rhs,
                'trK_rhs': trK_rhs,
                'Lambdabar_rhsU': Lambdabar_rhsU,
                'a_rhsDD': a_rhsDD,
                'alpha_rhs': alpha_rhs,
                'vet_rhsU': vet_rhsU,
                'bet_rhsU': bet_rhsU,
                'H': H,
                'MU': MU,
                'RbarDD': RbarDD
            }, {
                'h_rhsDD': Brhs.h_rhsDD,
                'cf_rhs': Brhs.cf_rhs,
                'trK_rhs': Brhs.trK_rhs,
                'Lambdabar_rhsU': Brhs.Lambdabar_rhsU,
                'a_rhsDD': Brhs.a_rhsDD,
                'alpha_rhs': gaugerhs.alpha_rhs,
                'vet_rhsU': gaugerhs.vet_rhsU,
                'bet_rhsU': gaugerhs.bet_rhsU,
                'H': bssncon.H,
                'MU': bssncon.MU,
                'RbarDD': Bq.RbarDD
            },
            suppress_message=True)
Ejemplo n.º 7
0
    def test_example_BSSN():
        import NRPy_param_funcs as par, reference_metric as rfm
        import BSSN.BSSN_RHSs as Brhs, BSSN.BSSN_quantities as Bq
        Parser.clear_namespace()
        parse(r"""
            % keydef basis [x, y, z]
            % ignore "\\%", "\qquad"
            % vardef 'deltaDD', 'vetU', 'lambdaU'
            % vardef -numeric -sym01 'hDD', 'aDD', 'RbarDD'
            % assign -numeric 'cf', 'alpha', 'trK', 'vetU', 'lambdaU'
            % parse \hat{\gamma}_{ij} = \delta_{ij}
            % assign -symbolic <H> -metric 'gammahatDD'
            % parse \bar{\gamma}_{ij} = h_{ij} + \hat{\gamma}_{ij}
            % assign -numeric -metric 'gammabarDD'
            \begin{align}
                %% replace LaTeX variable(s) with BSSN variable(s)
                % srepl "\bar{A}" -> "\text{a}", "\beta" -> "\text{vet}", "K" -> "\text{trK}", "\bar{\Lambda}" -> "\text{lambda}"
                % srepl "e^{{-4\phi}}" -> "\text{cf}^{{2}}"
                % srepl "\partial_t \phi = <1..> \\" -> "\text{cf_rhs} = -2 \text{cf} (<1..>) \\"
                % srepl -global "\partial_<1> \phi"         -> "\partial_<1> \text{cf} \frac{-1}{2 \text{cf}}"
                % srepl -global "\partial_<1> \text{phi}"   -> "\partial_<1> \text{cf} \frac{-1}{2 \text{cf}}"
                % srepl -global "\partial_<1> (\text{phi})" -> "\partial_<1> \text{cf} \frac{-1}{2 \text{cf}}"

                %% upwind pattern inside Lie derivative expansion
                % srepl -global "\text{vet}^<1> \partial_<1>" -> "\text{vet}^<1> \vphantom{upwind} \partial_<1>"

                %% enforce metric constraint gammabardet == gammahatdet
                % srepl -global "\bar{D}^i \bar{D}_k \text{vet}^k" -> "(\bar{D}^i \partial_k \text{vet}^k)"
                % srepl -global "\bar{D}_k \text{vet}^k" -> "(\partial_k \text{vet}^k + \frac{\partial_k \text{gammahatdet} \text{vet}^k}{2 \text{gammahatdet}})"

                % parse \bar{A}^i_j = \bar{\gamma}^{ik} \bar{A}_{kj}
                % srepl "\partial_t \bar{\gamma}" -> "\text{h_rhs}"
                \partial_t \bar{\gamma}_{ij} &= \mathcal{L}_\beta \bar{\gamma}_{ij}
                    + \frac{2}{3} \bar{\gamma}_{ij} \left(\alpha \bar{A}^k{}_k - \bar{D}_k \beta^k\right)
                    - 2 \alpha \bar{A}_{ij} \\

                \partial_t \phi &= \mathcal{L}_\beta \phi
                    + \frac{1}{6} \left(\bar{D}_k \beta^k - \alpha K \right) \\

                % parse \bar{A}^{ij} = \bar{\gamma}^{ik} \bar{\gamma}^{jl} \bar{A}_{kl}
                % srepl "\partial_t \text{trK}" -> "\text{trK_rhs}"
                \partial_t K &= \mathcal{L}_\beta K
                    + \frac{1}{3} \alpha K^{{2}}
                    + \alpha \bar{A}_{ij} \bar{A}^{ij}
                    - e^{{-4\phi}} \left(\bar{D}_i \bar{D}^i \alpha + 2 \bar{D}^i \alpha \bar{D}_i \phi\right) \\

                % parse \Pi^k_{ij} = \bar{\Gamma}^k_{ij} - \hat{\Gamma}^k_{ij}
                % parse \Pi_{ijk}  = \bar{\gamma}_{il} \Pi^l_{jk}
                % parse \Pi^k = \bar{\gamma}^{ij} \Pi^k_{ij}
                % srepl "\partial_t \text{lambda}" -> "\text{Lambdabar_rhs}"
                \partial_t \bar{\Lambda}^i &= \mathcal{L}_\beta \bar{\Lambda}^i + \bar{\gamma}^{jk} \hat{D}_j \hat{D}_k \beta^i
                    + \frac{2}{3} \Pi^i \bar{D}_k \beta^k + \frac{1}{3} \bar{D}^i \bar{D}_k \beta^k \\%
                    &\qquad- 2 \bar{A}^{ij} (\partial_j \alpha - 6 \alpha \partial_j \phi)
                    + 2 \alpha \bar{A}^{jk} \Pi^i_{jk} - \frac{4}{3} \alpha \bar{\gamma}^{ij} \partial_j K \\

                X_{ij} &= -2 \alpha \bar{D}_i \bar{D}_j \phi + 4 \alpha \bar{D}_i \phi \bar{D}_j \phi
                    + 2 \bar{D}_i \alpha \bar{D}_j \phi + 2 \bar{D}_j \alpha \bar{D}_i \phi
                    - \bar{D}_i \bar{D}_j \alpha + \alpha \bar{R}_{ij} \\
                \hat{X}_{ij} &= X_{ij} - \frac{1}{3} \bar{\gamma}_{ij} \bar{\gamma}^{kl} X_{kl} \\
                % srepl "\partial_t \text{a}" -> "\text{a_rhs}"
                \partial_t \bar{A}_{ij} &= \mathcal{L}_\beta \bar{A}_{ij}
                    - \frac{2}{3} \bar{A}_{ij} \bar{D}_k \beta^k
                    - 2 \alpha \bar{A}_{ik} \bar{A}^k_j
                    + \alpha \bar{A}_{ij} K
                    + e^{{-4\phi}} \hat{X}_{ij} \\

                \bar{R}_{ij} &= -\frac{1}{2} \bar{\gamma}^{kl} \hat{D}_k \hat{D}_l \bar{\gamma}_{ij}
                    + \frac{1}{2} (\bar{\gamma}_{ki} \hat{D}_j \bar{\Lambda}^k + \bar{\gamma}_{kj} \hat{D}_i \bar{\Lambda}^k)
                    + \frac{1}{2} \Pi^k (\Pi_{ijk} + \Pi_{jik}) \\%
                    &\qquad+ \bar{\gamma}^{kl} (\Pi^m_{ki} \Pi_{jml} + \Pi^m_{kj} \Pi_{iml} + \Pi^m_{ik} \Pi_{mjl}) 
            \end{align}
        """)
        par.set_parval_from_str('reference_metric::CoordSystem', 'Cartesian')
        par.set_parval_from_str('BSSN.BSSN_quantities::LeaveRicciSymbolic',
                                'True')
        rfm.reference_metric()
        Brhs.BSSN_RHSs()
        par.set_parval_from_str('BSSN.BSSN_quantities::LeaveRicciSymbolic',
                                'False')
        Bq.RicciBar__gammabarDD_dHatD__DGammaUDD__DGammaU()
        assert_equal(
            {
                'h_rhsDD': h_rhsDD,
                'cf_rhs': cf_rhs,
                'trK_rhs': trK_rhs,
                'Lambdabar_rhsU': Lambdabar_rhsU,
                'a_rhsDD': a_rhsDD,
                'RbarDD': RbarDD
            }, {
                'h_rhsDD': Brhs.h_rhsDD,
                'cf_rhs': Brhs.cf_rhs,
                'trK_rhs': Brhs.trK_rhs,
                'Lambdabar_rhsU': Brhs.Lambdabar_rhsU,
                'a_rhsDD': Brhs.a_rhsDD,
                'RbarDD': Bq.RbarDD
            })