#!/usr/bin/env python3 # # Aaron LI # Created: 2016-06-24 # Updated: 2016-07-04 # # Change logs: # 2016-07-04: # * Remove unnecessary configuration options # * Update radii unit to be "kpc", mass unit to be "Msun" # 2016-06-29: # * Update "plot()" to support electron number density profile # 2016-06-28: # * Implement plot function # * Adjust integration tolerances and progress report # * Fit smoothing splines to profiles using R `mgcv::gam()` # 2016-06-27: # * Implement potential profile calculation: # methods 'calc_density_total()' and 'calc_potential()' # 2016-06-26: # * Add document on gravitational potential calculation # 2016-06-25: # * Rename method 'interpolate()' to 'fit_spline()' # * Use 'InterpolatedUnivariateSpline' instead of 'interp1d' # * Implement 'calc_mass_total()' # * Update documentation # 2016-06-24: # * Update method 'gen_radius()' # """ Calculate the (gas and gravitational) mass profile and gravitational potential profile from the electron number density profile, under the assumption of hydrostatic equilibrium. The electron density profile and temperature profile are required. Assuming that the gas is in hydrostatic equilibrium with the gravitational potential and a spherically-symmetric distribution of the gas, we can write the hydrostatic equilibrium equation (HEE) of the ICM as (ref.[1], eq.(6)): derivative(P_gas, r) / rho_gas = - derivative(phi, r) = - G M_tot( R, i.e., outside the shell) The total gravitational potential may be considered to be the sum of the potentials of spherical shells of mass dM(r) = 4 * pi * rho(r) * r^2 dr, so we may calculate the gravitational potential at 'R' generated by an arbitrary spherically symmetric density distribution 'rho(r)' by adding the contributions to the potential produced by shells (1) with r < R, and (2) with r > R. In this way, we obtain phi(R) = - (G/R) * \int_0^R dM(r) - G * \int_R^{\inf} dM(r)/r = - 4*pi*G * ((1/R) * \int_0^R r^2 * rho(r) dr + \int_R^{\inf} r * rho(r) dr) ------------------------------------------------------------ References: [1] Ettori et al., 2013, Space Science Review, 177, 119-154 [2] Walker et al., 2012, MNRAS, 422, 3503 [3] Tremaine & Binney, Galactic Dynamics, 2nd edition, 2008 Sample configuration file: ------------------------------------------------------------ ## Configuration for `calc_mass_potential.py` ## Date: 2016-06-28 # electron density profile ne_profile = ne_profile.txt # cooling function profile cf_profile = coolfunc_profile.txt # temperature profile t_profile = t_profile.txt # number of data points for the output profiles (interpolation) num_dp = 1000 # output gas mass profile m_gas_profile = mass_gas_profile.txt m_gas_profile_image = mass_gas_profile.png # output total (gravitational) mass profile m_total_profile = mass_total_profile.txt m_total_profile_image = mass_total_profile.png # output total mass density profile rho_total_profile = rho_total_profile.txt rho_total_profile_image = rho_total_profile.png # output gravitational potential profile potential_profile = potential_profile.txt potential_profile_image = potential_profile.png ------------------------------------------------------------ """ import argparse import numpy as np import astropy.units as au import astropy.constants as ac import scipy.integrate as integrate from scipy.misc import derivative from configobj import ConfigObj import matplotlib.pyplot as plt from matplotlib.backends.backend_agg import FigureCanvasAgg as FigureCanvas from matplotlib.figure import Figure import rpy2.robjects as ro from rpy2.robjects.packages import importr from astro_params import AstroParams, ChandraPixel from projection import Projection plt.style.use("ggplot") class DensityProfile: """ Utilize the 3D (electron number or gas mass) density profile to calculate the following quantities: * 2D projected surface brightness (requires cooling function profile) * gas mass profile (integrated, M_gas( g/cm^3 c = au.solMass.to(au.g) / au.kpc.to(au.cm)**3 for i, r in enumerate(self.radius): dM_dr = derivative(lambda r: self.eval_spline("mass_total", r), r, dx=0.01) rho_total[i] = (dM_dr / (4 * np.pi * r**2)) * c self.rho_total = rho_total return rho_total def calc_potential(self, verbose=True): """ Calculate the gravitational potential profile. NOTE: The integral in the potential formula can only be integrated to the largest radius of availability. This limitation will affect the absolute values of the calculated potentials, but not the shape of the potential profile. """ def _int_inner(x): return x**2 * self.eval_spline(spline="rho_total", x=x) def _int_outer(x): return x * self.eval_spline(spline="rho_total", x=x) if self.rho_total is None: self.calc_density_total(verbose=verbose) if self.rho_total_spline is None: self.fit_spline(spline="rho_total", log10=True) potential = np.zeros(self.radius.shape) if verbose: print("Calculating the potential profile (#%d): ..." % len(self.radius), end="", flush=True) r_max = max(self.radius) c = - 4 * np.pi * au.kpc.to(au.cm)**2 for i, r in enumerate(self.radius): if verbose and (i+1) % 10 == 0: print("%d..." % (i+1), end="", flush=True) # total gravitational potential [ cm^2/s^2 ] potential[i] = c * ac.G.cgs.value * ( (1/r) * integrate.quad(_int_inner, a=0.0, b=r, epsrel=1e-2)[0] + integrate.quad(_int_outer, a=r, b=r_max, epsrel=1e-2)[0]) if verbose: print("DONE!", flush=True) self.potential = potential return potential def plot(self, profile, with_spline=True, ax=None, fig=None): x = self.radius xlabel = "3D Radius" xunit = "kpc" xscale = "log" yscale = "log" x_spl, y_spl = None, None if profile == "electron": x = self.r y = self.ne ylabel = "Electron number density" yunit = "cm$^{-3}$" if with_spline: x_spl = self.radius y_spl = self.eval_spline(spline="electron", x=self.radius) elif profile == "mass_gas": y = self.m_gas ylabel = "Gas mass" yunit = "M$_{\odot}$" elif profile == "mass_total": y = self.m_total ylabel = "Total mass" yunit = "M$_{\odot}$" elif profile == "rho_total": y = self.rho_total ylabel = "Total mass density" yunit = "g/cm$^3$" elif profile == "potential": y = self.potential ylabel = "Gravitational potential" yunit = "cm$^2$/s$^2$" yscale = "linear" else: raise ValueError("unknown profile name: %s" % profile) # if ax is None: fig, ax = plt.subplots(1, 1) ax.plot(x, y, linewidth=2) if with_spline and y_spl is not None: ax.plot(x_spl, y_spl, linewidth=2, linestyle="dashed") ax.set_xscale(xscale) ax.set_yscale(yscale) ax.set_xlim(min(x), max(x)) y_min, y_max = min(y), max(y) y_min = y_min/1.2 if y_min > 0 else y_min*1.2 y_max = y_max*1.2 if y_max > 0 else y_max/1.2 ax.set_ylim(y_min, y_max) ax.set_xlabel(r"%s (%s)" % (xlabel, xunit)) ax.set_ylabel(r"%s (%s)" % (ylabel, yunit)) fig.tight_layout() return (fig, ax) def save(self, profile, outfile): if profile == "mass_gas": data = np.column_stack([self.radius, self.radius_err, self.m_gas]) header = "radius[kpc] radius_err[kpc] mass_gas(