From 463189a4d4bfe47ab12b9b8f798a876e97cf97b5 Mon Sep 17 00:00:00 2001 From: knarf Date: Mon, 1 Mar 2010 22:07:39 +0000 Subject: compute baryonic mass and proper distance, simply code a bit git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/TOVSolver/trunk@81 1bdb13ef-5d69-4035-bb54-08abeb3aa7f1 --- src/tov.c | 87 +++++++++++++++++++++++++++-------------------------------- src/utils.inc | 8 ++++++ 2 files changed, 48 insertions(+), 47 deletions(-) diff --git a/src/tov.c b/src/tov.c index f136c38..5aedb0d 100644 --- a/src/tov.c +++ b/src/tov.c @@ -17,6 +17,8 @@ #include "tov.h" +#define NUMVARS 6 + #define velx (&vel[0*cctk_lsh[0]*cctk_lsh[1]*cctk_lsh[2]]) #define vely (&vel[1*cctk_lsh[0]*cctk_lsh[1]*cctk_lsh[2]]) #define velz (&vel[2*cctk_lsh[0]*cctk_lsh[1]*cctk_lsh[2]]) @@ -46,6 +48,8 @@ CCTK_REAL * TOV_rbar_1d=0; CCTK_REAL * TOV_press_1d=0; CCTK_REAL * TOV_phi_1d=0; CCTK_REAL * TOV_m_1d=0; +CCTK_REAL * TOV_mbary_1d=0; +CCTK_REAL * TOV_rprop_1d=0; #include "utils.inc" @@ -115,10 +119,10 @@ void TOV_C_ParamCheck(CCTK_ARGUMENTS) @endhistory @@*/ void TOV_C_Source_RHS(CCTK_REAL r, CCTK_REAL K, CCTK_REAL Gamma, - CCTK_REAL old_data[4], CCTK_REAL source_data[4]) + CCTK_REAL old_data[NUMVARS], CCTK_REAL source_data[NUMVARS]) { CCTK_REAL LOCAL_TINY, PI; - CCTK_REAL press, rho, eps, mu, m, phi; + CCTK_REAL press, rho, eps, mu, m, phi, mbary, rprop; CCTK_REAL log_rbar_over_r, r_minus_two_m; LOCAL_TINY = 1.0e-35; @@ -130,6 +134,8 @@ void TOV_C_Source_RHS(CCTK_REAL r, CCTK_REAL K, CCTK_REAL Gamma, m = old_data[1]; phi = old_data[2]; log_rbar_over_r = old_data[3]; + mbary = old_data[4]; + rprop = old_data[5]; rho = pow(press / K, 1.0 / Gamma); eps = press / (Gamma - 1.0) / rho; @@ -139,10 +145,11 @@ void TOV_C_Source_RHS(CCTK_REAL r, CCTK_REAL K, CCTK_REAL Gamma, if ((r==0.0) && (m==0.0)) { - /* source_data[0] = 0.0; */ source_data[1] = 0.0; source_data[2] = 0.0; source_data[3] = 0.0; + source_data[4] = 0.0; + source_data[5] = 0.0; } else { @@ -152,6 +159,8 @@ void TOV_C_Source_RHS(CCTK_REAL r, CCTK_REAL K, CCTK_REAL Gamma, (m + 4*PI * r*r*r * press) / r_minus_two_m / r; source_data[1] = 4*PI * r*r * mu; source_data[3] = (sqrt(r) - sqrt(r_minus_two_m)) / r / sqrt(r_minus_two_m); + source_data[5] = 1.0/sqrt(1.0-2.0*m/r); + source_data[4] = source_data[5] * 4*PI * rho * r*r; } } @@ -168,7 +177,6 @@ void TOV_C_Source_RHS(CCTK_REAL r, CCTK_REAL K, CCTK_REAL Gamma, @history @endhistory @@*/ - void TOV_C_Integrate_RHS(CCTK_ARGUMENTS) { DECLARE_CCTK_ARGUMENTS @@ -177,8 +185,9 @@ void TOV_C_Integrate_RHS(CCTK_ARGUMENTS) CCTK_REAL LOCAL_TINY; CCTK_INT star, star_i, i, TOV_Surface_Index; - CCTK_REAL old_data[4], source_data[4], in_data[4], new_data[4], - k1[4], k2[4], k3[4], k4[4]; + CCTK_REAL old_data[NUMVARS], source_data[NUMVARS], + in_data[NUMVARS], new_data[NUMVARS], + k1[NUMVARS], k2[NUMVARS], k3[NUMVARS], k4[NUMVARS]; CCTK_REAL Surface_Mass, factor, local_rho; CCTK_REAL rho_central; @@ -193,6 +202,8 @@ void TOV_C_Integrate_RHS(CCTK_ARGUMENTS) assert(TOV_press_1d!=0); assert(TOV_phi_1d!=0); assert(TOV_m_1d!=0); + assert(TOV_mbary_1d!=0); + assert(TOV_rprop_1d!=0); /* do it for all stars */ for (star=0; star < TOV_Num_TOVs; star++) @@ -240,6 +251,8 @@ void TOV_C_Integrate_RHS(CCTK_ARGUMENTS) TOV_C_fill(&(TOV_phi_1d [star_i]), TOV_Num_Radial, 0.0); TOV_C_fill(&(TOV_rbar_1d [star_i]), TOV_Num_Radial, 0.0); TOV_C_fill(&(TOV_r_1d [star_i]), TOV_Num_Radial, 0.0); + TOV_C_fill(&(TOV_mbary_1d[star_i]), TOV_Num_Radial, 0.0); + TOV_C_fill(&(TOV_rprop_1d[star_i]), TOV_Num_Radial, 0.0); /* set start values */ TOV_press_1d[star_i] = TOV_K[star] * @@ -254,10 +267,12 @@ void TOV_C_Integrate_RHS(CCTK_ARGUMENTS) TOV_Surface[star] = -1.0; TOV_Surface_Index = -1.0; +#define RKLOOP for (int rk=0; rk