def L3_integrate( num_prof, R, M, P_s, T_s, rho_s, R1, R2, mat_id_L1, T_rho_type_id_L1, T_rho_args_L1, mat_id_L2, T_rho_type_id_L2, T_rho_args_L2, mat_id_L3, T_rho_type_id_L3, T_rho_args_L3, ): """Integration of a 2 layer spherical planet. Parameters ---------- num_prof : int Number of profile integration steps. R : float Radii of the planet (m). M : float Mass of the planet (kg). P_s : float Pressure at the surface (Pa). T_s : float Temperature at the surface (K). rho_s : float Density at the surface (kg m^-3). R1 : float Boundary between layers 1 and 2 (m). R2 : float Boundary between layers 2 and 3 (m). mat_id_L1 : int Material id for layer 1. T_rho_type_id_L1 : int Relation between A1_T and A1_rho to be used in layer 1. T_rho_args_L1 : [float] Extra arguments to determine the relation in layer 1. mat_id_L2 : int Material id for layer 2. T_rho_type_id_L2 : int Relation between A1_T and A1_rho to be used in layer 2. T_rho_args_L2 : [float] Extra arguments to determine the relation in layer 2. mat_id_L3 : int Material id for layer 3. T_rho_type_id_L3 : int Relation between A1_T and A1_rho to be used in layer 3. T_rho_args_L3 : [float] Extra arguments to determine the relation in layer 3. Returns ------- A1_r : [float] The profile radii, in increasing order (m). A1_m_enc : [float] The cummulative mass at each profile radius (kg). A1_P : [float] The pressure at each profile radius (Pa). A1_T : [float] The temperature at each profile radius (K). A1_rho : [float] The density at each profile radius (kg m^-3). A1_u : [float] The specific internal energy at each profile radius (J kg^-1). A1_mat_id : [float] The ID of the material at each profile radius. """ A1_r = np.linspace(R, 0, int(num_prof)) A1_m_enc = np.zeros(A1_r.shape) A1_P = np.zeros(A1_r.shape) A1_T = np.zeros(A1_r.shape) A1_rho = np.zeros(A1_r.shape) A1_u = np.zeros(A1_r.shape) A1_mat_id = np.zeros(A1_r.shape) u_s = eos.u_rho_T(rho_s, T_s, mat_id_L3) T_rho_args_L3 = set_T_rho_args(T_s, rho_s, T_rho_type_id_L3, T_rho_args_L3, mat_id_L3) dr = A1_r[0] - A1_r[1] A1_m_enc[0] = M A1_P[0] = P_s A1_T[0] = T_s A1_rho[0] = rho_s A1_u[0] = u_s A1_mat_id[0] = mat_id_L3 for i in range(1, A1_r.shape[0]): # Layer 3 if A1_r[i] > R2: rho = A1_rho[i - 1] mat_id = mat_id_L3 T_rho_type_id = T_rho_type_id_L3 T_rho_args = T_rho_args_L3 rho0 = rho # Layer 2, 3 boundary elif A1_r[i] <= R2 and A1_r[i - 1] > R2: # New density, continuous temperature unless fixed entropy if T_rho_type_id_L2 == gv.type_ent: rho = eos.find_rho( A1_P[i - 1], mat_id_L2, T_rho_type_id_L2, T_rho_args_L2, A1_rho[i - 1], 1e5, ) else: rho = eos.rho_P_T(A1_P[i - 1], A1_T[i - 1], mat_id_L2) T_rho_args_L2 = set_T_rho_args(A1_T[i - 1], rho, T_rho_type_id_L2, T_rho_args_L2, mat_id_L2) mat_id = mat_id_L2 T_rho_type_id = T_rho_type_id_L2 T_rho_args = T_rho_args_L2 rho0 = A1_rho[i - 1] # Layer 2 elif A1_r[i] > R1: rho = A1_rho[i - 1] mat_id = mat_id_L2 T_rho_type_id = T_rho_type_id_L2 T_rho_args = T_rho_args_L2 rho0 = rho # Layer 1, 2 boundary elif A1_r[i] <= R1 and A1_r[i - 1] > R1: # New density, continuous temperature unless fixed entropy if T_rho_type_id_L1 == gv.type_ent: rho = eos.find_rho( A1_P[i - 1], mat_id_L1, T_rho_type_id_L1, T_rho_args_L1, A1_rho[i - 1], 1e5, ) else: rho = eos.rho_P_T(A1_P[i - 1], A1_T[i - 1], mat_id_L1) T_rho_args_L1 = set_T_rho_args(A1_T[i - 1], rho, T_rho_type_id_L1, T_rho_args_L1, mat_id_L1) mat_id = mat_id_L1 T_rho_type_id = T_rho_type_id_L1 T_rho_args = T_rho_args_L1 rho0 = A1_rho[i - 1] # Layer 1 elif A1_r[i] <= R1: rho = A1_rho[i - 1] mat_id = mat_id_L1 T_rho_type_id = T_rho_type_id_L1 T_rho_args = T_rho_args_L1 rho0 = A1_rho[i - 1] A1_m_enc[i] = A1_m_enc[i - 1] - 4 * np.pi * A1_r[i - 1]**2 * rho * dr A1_P[i] = A1_P[i - 1] + gv.G * A1_m_enc[i - 1] * rho / (A1_r[i - 1]** 2) * dr A1_rho[i] = eos.find_rho(A1_P[i], mat_id, T_rho_type_id, T_rho_args, rho0, 1.1 * rho) A1_T[i] = T_rho(A1_rho[i], T_rho_type_id, T_rho_args, mat_id) A1_u[i] = eos.u_rho_T(A1_rho[i], A1_T[i], mat_id) A1_mat_id[i] = mat_id # Update the T-rho parameters if mat_id == gv.id_HM80_HHe and T_rho_type_id == gv.type_adb: T_rho_args = set_T_rho_args(A1_T[i], A1_rho[i], T_rho_type_id, T_rho_args, mat_id) if A1_m_enc[i] < 0: break return A1_r, A1_m_enc, A1_P, A1_T, A1_rho, A1_u, A1_mat_id
def L1_integrate_out( r, dr, m_enc, P, T, u, mat_id, T_rho_type_id, T_rho_args, rho_min=0, P_min=0 ): """Integrate a new layer of a spherical planet outwards. Parameters ---------- r : float The radius at the base (m). dr : float The radius step for the integration (m). m_enc : float The enclosed mass at the base (Pa). P : float The pressure at the base (Pa). T : float The temperature at the base (K). u : float The specific internal energy at the base (J kg^-1). mat_id : int The material ID. T_rho_type_id : int The ID for the temperature-density relation. T_rho_args : [float] Extra arguments for the temperature-density relation. rho_min : float The minimum density (must be >= 0) at which the new layer will stop. P_min : float The minimum pressure (must be >= 0) at which the new layer will stop. Returns ------- A1_r : [float] The profile radii, in increasing order (m). A1_m_enc : [float] The cummulative mass at each profile radius (kg). A1_P : [float] The pressure at each profile radius (Pa). A1_T : [float] The temperature at each profile radius (K). A1_rho : [float] The density at each profile radius (kg m^-3). A1_u : [float] The specific internal energy at each profile radius (J kg^-1). A1_mat_id : [float] The ID of the material at each profile radius. """ # Initialise the profile arrays A1_r = [r] A1_m_enc = [m_enc] A1_P = [P] A1_T = [T] A1_u = [u] A1_mat_id = [mat_id] A1_rho = [eos.rho_P_T(A1_P[0], A1_T[0], mat_id)] # Integrate outwards until the minimum density (or zero pressure) while A1_rho[-1] > rho_min and A1_P[-1] > P_min: A1_r.append(A1_r[-1] + dr) A1_m_enc.append( A1_m_enc[-1] + 4 * np.pi * A1_r[-1] * A1_r[-1] * A1_rho[-1] * dr ) A1_P.append(A1_P[-1] - gv.G * A1_m_enc[-1] * A1_rho[-1] / (A1_r[-1] ** 2) * dr) if A1_P[-1] <= 0: # Add dummy values which will be removed along with the -ve P A1_rho.append(0) A1_T.append(0) A1_u.append(0) A1_mat_id.append(0) break # Update the T-rho parameters if T_rho_type_id == gv.type_adb and mat_id == gv.id_HM80_HHe: T_rho_args = set_T_rho_args( A1_T[-1], A1_rho[-1], T_rho_type_id, T_rho_args, mat_id, ) rho = eos.find_rho( A1_P[-1], mat_id, T_rho_type_id, T_rho_args, 0.9 * A1_rho[-1], A1_rho[-1], ) A1_rho.append(rho) A1_T.append(T_rho(rho, T_rho_type_id, T_rho_args, mat_id)) A1_u.append(eos.u_rho_T(rho, A1_T[-1], mat_id)) A1_mat_id.append(mat_id) # Remove the duplicate first step and the final too-low density or too-low # pressure step return ( A1_r[1:-1], A1_m_enc[1:-1], A1_P[1:-1], A1_T[1:-1], A1_rho[1:-1], A1_u[1:-1], A1_mat_id[1:-1], )
def L1_integrate(num_prof, R, M, P_s, T_s, rho_s, mat_id, T_rho_type_id, T_rho_args): """Integration of a 1 layer spherical planet. Parameters ---------- num_prof : int Number of profile integration steps. R : float Radii of the planet (m). M : float Mass of the planet (kg). P_s : float Pressure at the surface (Pa). T_s : float Temperature at the surface (K). rho_s : float Density at the surface (kg m^-3). mat_id : int Material id. T_rho_type_id : int Relation between A1_T and A1_rho to be used. T_rho_args : [float] Extra arguments to determine the relation. Returns ------- A1_r : [float] The profile radii, in increasing order (m). A1_m_enc : [float] The cummulative mass at each profile radius (kg). A1_P : [float] The pressure at each profile radius (Pa). A1_T : [float] The temperature at each profile radius (K). A1_rho : [float] The density at each profile radius (kg m^-3). A1_u : [float] The specific internal energy at each profile radius (J kg^-1). A1_mat_id : [float] The ID of the material at each profile radius. """ # Initialise the profile arrays A1_r = np.linspace(R, 0, int(num_prof)) A1_m_enc = np.zeros(A1_r.shape) A1_P = np.zeros(A1_r.shape) A1_T = np.zeros(A1_r.shape) A1_rho = np.zeros(A1_r.shape) A1_u = np.zeros(A1_r.shape) A1_mat_id = np.ones(A1_r.shape) * mat_id u_s = eos.u_rho_T(rho_s, T_s, mat_id) # Set the T-rho relation parameters T_rho_args = set_T_rho_args(T_s, rho_s, T_rho_type_id, T_rho_args, mat_id) dr = A1_r[0] - A1_r[1] # Set the surface values A1_m_enc[0] = M A1_P[0] = P_s A1_T[0] = T_s A1_rho[0] = rho_s A1_u[0] = u_s # Integrate inwards for i in range(1, A1_r.shape[0]): A1_m_enc[i] = ( A1_m_enc[i - 1] - 4 * np.pi * A1_r[i - 1] ** 2 * A1_rho[i - 1] * dr ) A1_P[i] = ( A1_P[i - 1] + gv.G * A1_m_enc[i - 1] * A1_rho[i - 1] / (A1_r[i - 1] ** 2) * dr ) A1_rho[i] = eos.find_rho( A1_P[i], mat_id, T_rho_type_id, T_rho_args, A1_rho[i - 1], 1.1 * A1_rho[i - 1], ) A1_T[i] = T_rho(A1_rho[i], T_rho_type_id, T_rho_args, mat_id) A1_u[i] = eos.u_rho_T(A1_rho[i], A1_T[i], mat_id) # Update the T-rho parameters if mat_id == gv.id_HM80_HHe and T_rho_type_id == gv.type_adb: T_rho_args = set_T_rho_args( A1_T[i], A1_rho[i], T_rho_type_id, T_rho_args, mat_id ) # Stop if run out of mass if A1_m_enc[i] < 0: return A1_r, A1_m_enc, A1_P, A1_T, A1_rho, A1_u, A1_mat_id return A1_r, A1_m_enc, A1_P, A1_T, A1_rho, A1_u, A1_mat_id
def L1_rho_eq_po_from_V( A1_r_eq, A1_V_eq, A1_r_po, A1_V_po, P_0, P_s, rho_0, rho_s, mat_id_L1, T_rho_type_id_L1, T_rho_args_L1, ): """Compute densities of equatorial and polar profiles given the potential for a 1 layer planet. Parameters ---------- A1_r_eq : [float] Points at equatorial profile where the solution is defined (m). A1_V_eq : [float] Equatorial profile of potential (J). A1_r_po : [float] Points at equatorial profile where the solution is defined (m). A1_V_po : [float] Polar profile of potential (J). P_0 : float Pressure at the center of the planet (Pa). P_s : float Pressure at the surface of the planet (Pa). rho_0 : float Density at the center of the planet (kg m^-3). rho_s : float Density at the surface of the planet (kg m^-3). mat_id_L1 : int Material id for layer 1. T_rho_type_id_L1 : int Relation between T and rho to be used in layer 1. T_rho_args_L1 : [float] Extra arguments to determine the relation in layer 1. Returns ------- A1_rho_eq : [float] Equatorial profile of densities (kg m^-3). A1_rho_po : [float] Polar profile of densities (kg m^-3). """ A1_P_eq = np.zeros(A1_V_eq.shape[0]) A1_P_po = np.zeros(A1_V_po.shape[0]) A1_rho_eq = np.zeros(A1_V_eq.shape[0]) A1_rho_po = np.zeros(A1_V_po.shape[0]) A1_P_eq[0] = P_0 A1_P_po[0] = P_0 A1_rho_eq[0] = rho_0 A1_rho_po[0] = rho_0 # equatorial profile for i in range(A1_r_eq.shape[0] - 1): gradV = A1_V_eq[i + 1] - A1_V_eq[i] gradP = -A1_rho_eq[i] * gradV A1_P_eq[i + 1] = A1_P_eq[i] + gradP # avoid overspin if A1_P_eq[i + 1] > A1_P_eq[i]: A1_rho_eq[i + 1:] = A1_rho_eq[i] break # compute density if A1_P_eq[i + 1] >= P_s: A1_rho_eq[i + 1] = eos.find_rho( A1_P_eq[i + 1], mat_id_L1, T_rho_type_id_L1, T_rho_args_L1, rho_s * 0.1, A1_rho_eq[i], ) else: A1_rho_eq[i + 1] = 0.0 break # polar profile for i in range(A1_r_po.shape[0] - 1): gradV = A1_V_po[i + 1] - A1_V_po[i] gradP = -A1_rho_po[i] * gradV A1_P_po[i + 1] = A1_P_po[i] + gradP # avoid overspin if A1_P_po[i + 1] > A1_P_po[i]: A1_rho_po[i + 1:] = A1_rho_po[i] break # compute density if A1_P_po[i + 1] >= P_s: A1_rho_po[i + 1] = eos.find_rho( A1_P_po[i + 1], mat_id_L1, T_rho_type_id_L1, T_rho_args_L1, rho_s * 0.1, A1_rho_po[i], ) else: A1_rho_po[i + 1] = 0.0 break return A1_rho_eq, A1_rho_po