#!/usr/bin/env python3 # # Weitian LI # Created: 2016-06-24 # Updated: 2016-06-27 # # Change logs: # 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-24 # redshift used for pixel to distance conversion redshift = # electron density profile ne_profile = ne_profile.txt # cooling function profile cf_profile = coolfunc_profile.txt # unit of the CF profile radius (default: pixel) cf_unit = "pixel" # temperature profile t_profile = t_profile.txt # unit of the T profile radius (default: pixel) t_unit = "pixel" # number of data points for the output profile calculation num_dp = 1000 # output gas mass profile m_gas_profile = mass_gas_profile.txt # output total (gravitational) mass profile m_total_profile = mass_total_profile.txt # output total mass density profile rho_total_profile = rho_total_profile.txt # output gravitational potential profile potential_profile = potential_profile.txt ------------------------------------------------------------ """ import argparse import numpy as np import astropy.units as au import astropy.constants as ac import scipy.interpolate as interpolate import scipy.integrate as integrate from scipy.misc import derivative from configobj import ConfigObj from astro_params import AstroParams, ChandraPixel from projection import Projection 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(>> Fitting spline to total mass profile ...") self.m_total_spline = interpolate.InterpolatedUnivariateSpline( x=self.radius, y=self.m_total, ext="const") for i, r in enumerate(self.radius): rho_total[i] = (derivative(self.m_total_spline, r, dx=0.01*au.kpc.to(au.cm)) / (4 * np.pi * r**2)) self.rho_total = rho_total if verbose: print(">>> Fitting spline to total mass density profile ...") self.rho_total_spline = interpolate.InterpolatedUnivariateSpline( x=self.radius, y=rho_total, ext="const") 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.rho_total_spline(x) def _int_outer(x): return x * self.rho_total_spline(x) if self.rho_total is None: self.calc_density_total(verbose=verbose) 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) for i, r in enumerate(self.radius): if verbose and (i+1) % 50 == 0: print("%d..." % (i+1), end="", flush=True) potential[i] = - 4 * np.pi * ac.G.cgs.value * ( (1/r) * integrate.quad(_int_inner, a=0.0, b=r, epsabs=1e-5*au.kpc.to(au.cm), epsrel=1e-3)[0] + integrate.quad(_int_outer, a=r, b=r_max, epsabs=1e-5*au.kpc.to(au.cm), epsrel=1e-3)[0]) if verbose: print("DONE!", flush=True) self.potential = potential return potential def plot(self, profile, ax=None, fig=None): pass def save(self, profile, outfile): if profile == "mass_gas": data = np.column_stack([self.radius, self.radius_err, self.m_gas]) header = "radius[cm] radius_err[cm] mass_gas(