def update(self, potential, af): # the main procedure for recomputing the density generated by DF # and representing it by either a spherical- or azimuthal-harmonic expansion if not hasattr(self, 'df'): return # can only update DF-based components gm = agama.GalaxyModel(potential, self.df, af) densfnc = lambda x: gm.moments(x, dens=True, vel=False, vel2=False) if self.disklike: self.density = agama.Density(density=densfnc, type='DensityAzimuthalHarmonic', gridsizer=self.sizeradialcyl, rmin=self.rmincyl, rmax=self.rmaxcyl, gridsizez=self.sizeverticalcyl, zmin=self.zmincyl, zmax=self.zmaxcyl, mmax=0, symmetry='a') else: self.density = agama.Density(density=densfnc, type='DensitySphericalHarmonic', gridsizer=self.sizeradialsph, rmin=self.rminsph, rmax=self.rmaxsph, lmax=self.lmaxangularsph, mmax=0, symmetry='a')
def __init__(self, filename): ''' Initialize starting values of scaled parameters and and define upper/lower limits on them; also obtain the true values by analyzing the input file name ''' self.initValues = numpy.array( [9., 0., 0.5, 4.0, 1.0, 6.0, 1.5, 1.0, 1.0, 1.0, 1.]) self.minValues = numpy.array( [6., -1., 0.0, 2.2, .25, 3.2, 0.0, 0.5, 0.1, 0.1, -1.]) self.maxValues = numpy.array( [10., 2., 1.9, 6.0, 2.5, 12., 2.8, 2.0, 2.8, 2.8, 3.]) self.labels = ( \ r'$\log_{10}(\rho_1)$', r'$\log_{10}(R_{scale})$', r'$\gamma$', r'$\beta$', r'$\alpha$', \ r'$B_{DF}$', r'$\Gamma_{DF}$', r'$\eta$', r'$g_r$', r'$h_r$', r'$\log_{10}(J_0)$') self.numPotParams = 5 # potential params come first in the list self.numDFParams = 6 # DF params come last in the list self.scaleRadius = 1. # fixed radius r_1 at which the DM density is constrained (rho_1) # parse true values m = re.match( r'gs(\d+)_bs(\d+)_rcrs([a-z\d]+)_rarc([a-z\d]+)_([a-z]+)_(\d+)mpc3', filename) n = re.match(r'data_([ch])_rh(\d+)_rs(\d+)_gs(\d+)', filename) if m: self.truePotential = agama.Potential( type='SpheroidDensity', densityNorm=float(m.group(6)) * 1e6, scaleRadius=1.0, gamma=1.0 if m.group(5) == 'cusp' else 0.0, beta=3.0) self.scaleRadius = float(m.group(3)) * 0.01 self.tracerParams = dict(type='SpheroidDensity', densityNorm=1.0, scaleRadius=self.scaleRadius, gamma=float(m.group(1)) * 0.01, beta=float(m.group(2)) * 0.1) # normalize the tracer density profile to have a total mass of unity self.tracerParams["densityNorm"] /= agama.Density( **self.tracerParams).totalMass() self.tracerDensity = agama.Density(**self.tracerParams) elif n: self.truePotential = agama.Potential( type='SpheroidDensity', densityNorm=3.021516e7 if n.group(1) == 'c' else 2.387329e7, scaleRadius=float(n.group(2)), gamma=0.0 if n.group(1) == 'c' else 1.0, beta=4.0) self.scaleRadius = float(n.group(3)) * 0.01 self.tracerParams = dict(type='SpheroidDensity', densityNorm=1.0, scaleRadius=self.scaleRadius, gamma=float(m.group(4)) * 0.1, beta=5.0) self.tracerParams["densityNorm"] /= agama.Density( **self.tracerParams).totalMass() self.tracerDensity = agama.Density(**self.tracerParams) else: print "Can't determine true parameters!"
def createModel(self, params): ''' create a model (potential and DF) specified by the given [scaled] parameters ''' potential = agama.Potential(type='Spheroid', densityNorm=10**params[0], scaleRadius=10**params[1], gamma=params[2], beta=params[3], alpha=params[4]) density = agama.Density(type='Spheroid', scaleRadius=10**params[self.numPotParams + 0], gamma=params[self.numPotParams + 1], beta=params[self.numPotParams + 2], alpha=params[self.numPotParams + 3]) df = agama.DistributionFunction(type='QuasiSpherical', potential=potential, density=density, beta0=params[self.numPotParams + 4], r_a=10**params[self.numPotParams + 5]) # check if the DF is everywhere nonnegative j = numpy.logspace(-5, 10, 200) if any(df(numpy.column_stack((j, j * 0 + 1e-10, j * 0))) <= 0): raise ValueError("Bad DF") return potential, df
def makeXBar(**params): densityNorm = params['densityNorm'] x0 = params['x0'] y0 = params['y0'] z0 = params['z0'] xc = params['xc'] yc = params['yc'] c = params['c'] alpha = params['alpha'] cpar = params['cpar'] cperp = params['cperp'] m = params['m'] n = params['n'] outerCutoffRadius = params['outerCutoffRadius'] def density(xyz): r = numpy.sum(xyz**2, axis=1)**0.5 a = (((abs(xyz[:, 0]) / x0)**cperp + (abs(xyz[:, 1]) / y0)**cperp)**(cpar / cperp) + (abs(xyz[:, 2]) / z0)**cpar)**(1 / cpar) ap = (((xyz[:, 0] + c * xyz[:, 2]) / xc)**2 + (xyz[:, 1] / yc)**2)**(0.5) am = (((xyz[:, 0] - c * xyz[:, 2]) / xc)**2 + (xyz[:, 1] / yc)**2)**(0.5) return (densityNorm / numpy.cosh(a**m) * numpy.exp(-(r / outerCutoffRadius)**2) * (1 + alpha * (numpy.exp(-ap**n) + numpy.exp(-am**n)))) return agama.Density(density)
def makePotentialModel(params): # combined 4 components and the CMC represented by a single triaxial CylSpline potential mmax = 12 # order of azimuthal Fourier expansion (higher order means better accuracy, # but values greater than 12 *significantly* slow down the computation!) pot_bary = agama.Potential(type='CylSpline', density=agama.Density( makeDensityModel(params), makeCMC(0.2e10, 0.25, 0.05, 0.5)), symmetry='t', mmax=mmax, gridsizeR=25, gridsizez=25, Rmin=0.1, Rmax=40, zmin=0.05, zmax=20) # flattened axisymmetric dark halo with the Einasto profile pot_dark = agama.Potential(type='Multipole', density='Spheroid', axisratioz=0.8, gamma=0, beta=0, outerCutoffRadius=1.84, cutoffStrength=0.74, densityNorm=0.0263e10, gridsizer=26, rmin=0.01, rmax=1000, lmax=8) return agama.Potential(pot_bary, pot_dark)
def updatePotential(self): print("Updating potential...") # sort out density and potential components into several groups densitySph = [] densityCyl = [] potentials = [] for component in self.components: dens = component.getDensity() if dens is not None: if component.disklike: densityCyl.append(dens) else: densitySph.append(dens) else: potentials.append(component.getPotential()) # create a single Multipole potential for all non-disk-like density components if len(densitySph) > 0: potentials.append( agama.Potential(type='Multipole', density=agama.Density(*densitySph), gridsizer=self.sizeradialsph, rmin=self.rminsph, rmax=self.rmaxsph, lmax=self.lmaxangularsph, mmax=0, symmetry='a')) # create a single CylSpline potential representing all disk-like density components if len(densityCyl) > 0: potentials.append( agama.Potential(type='CylSpline', density=agama.Density(*densityCyl), gridsizer=self.sizeradialcyl, rmin=self.rmincyl, rmax=self.rmaxcyl, gridsizez=self.sizeverticalcyl, zmin=self.zmincyl, zmax=self.zmaxcyl, mmax=0, symmetry='a')) # combine all potential components and reinitialize the action finder self.potential = agama.Potential(*potentials) print("Updating action finder...") self.af = agama.ActionFinder(self.potential)
def makeDensityModel(params): ind = 0 densityDisk = makeDisk(surfaceDensity=params[ind + 0], scaleRadius=params[ind + 1], innerCutoffRadius=params[ind + 2], scaleHeight=params[ind + 3], sersicIndex=params[ind + 4], verticalSersicIndex=params[ind + 5]) ind += 6 densityXBar = makeXBar(densityNorm=params[ind + 0], x0=params[ind + 1], y0=params[ind + 2], z0=params[ind + 3], cpar=params[ind + 4], cperp=params[ind + 5], m=params[ind + 6], outerCutoffRadius=params[ind + 7], alpha=params[ind + 8], c=params[ind + 9], n=params[ind + 10], xc=params[ind + 11], yc=params[ind + 12]) ind += 13 densityLongBar1 = makeLongBar(densityNorm=params[ind + 0], x0=params[ind + 1], y0=params[ind + 2], scaleHeight=params[ind + 3], cperp=params[ind + 4], cpar=params[ind + 5], outerCutoffRadius=params[ind + 6], innerCutoffRadius=params[ind + 7], outerCutoffStrength=params[ind + 8], innerCutoffStrength=params[ind + 9]) ind += 10 densityLongBar2 = makeLongBar(densityNorm=params[ind + 0], x0=params[ind + 1], y0=params[ind + 2], scaleHeight=params[ind + 3], cperp=params[ind + 4], cpar=params[ind + 5], outerCutoffRadius=params[ind + 6], innerCutoffRadius=params[ind + 7], outerCutoffStrength=params[ind + 8], innerCutoffStrength=params[ind + 9]) ind += 10 assert len(params) == ind, 'invalid number of parameters' return agama.Density(densityDisk, densityXBar, densityLongBar1, densityLongBar2)
def makeDisk(**params): surfaceDensity = params['surfaceDensity'] scaleRadius = params['scaleRadius'] scaleHeight = params['scaleHeight'] innerCutoffRadius = params['innerCutoffRadius'] sersicIndex = params['sersicIndex'] verticalSersicIndex = params['verticalSersicIndex'] def density(xyz): R = (xyz[:, 0]**2 + xyz[:, 1]**2)**0.5 return (surfaceDensity / (4 * scaleHeight) * numpy.exp(-(R / scaleRadius)**sersicIndex - innerCutoffRadius / (R + 1e-100)) / numpy.cosh( (abs(xyz[:, 2]) / scaleHeight)**verticalSersicIndex)) return agama.Density(density)
def makeLongBar(**params): densityNorm = params['densityNorm'] x0 = params['x0'] y0 = params['y0'] cpar = params['cpar'] cperp = params['cperp'] scaleHeight = params['scaleHeight'] innerCutoffRadius = params['innerCutoffRadius'] outerCutoffRadius = params['outerCutoffRadius'] innerCutoffStrength = params['innerCutoffStrength'] outerCutoffStrength = params['outerCutoffStrength'] def density(xyz): R = (xyz[:, 0]**2 + xyz[:, 1]**2)**0.5 a = ((abs(xyz[:, 0]) / x0)**cperp + (abs(xyz[:, 1]) / y0)**cperp)**(1 / cperp) return densityNorm / numpy.cosh( xyz[:, 2] / scaleHeight)**2 * numpy.exp( -a**cpar - (R / outerCutoffRadius)**outerCutoffStrength - (innerCutoffRadius / R)**innerCutoffStrength) return agama.Density(density)
#!/usr/bin/python """ This example illustrates the use of Schwarzschild method to construct a triaxial Dehnen model. """ import agama, numpy # density profile den = agama.Density(type='Dehnen', axisRatioY=0.8, axisRatioZ=0.6, scaleRadius=1, mass=1, gamma=1) # potential corresponding to this density pot = agama.Potential(type='Multipole', lmax=8, mmax=6, density=den) # target1 is discretized density profile target1 = agama.Target(type='DensityClassicLinear', gridr=agama.nonuniformGrid(25, 0.1, 20.0), \ axisratioy=0.8, axisratioz=0.6, stripsPerPane=2) # target2 is discretized kinematic constraint (we will enforce velocity isotropy) target2 = agama.Target(type='KinemShell', gridR=agama.nonuniformGrid(15, 0.2, 10.0), degree=1) # construct initial conditions for the orbit library initcond, weightprior = den.sample(5000, potential=pot) # integration time is 50 orbital periods inttimes = 50 * pot.Tcirc(initcond)
def isFloat(obj): return isinstance(obj, float) def isTuple(obj, length): return isinstance(obj, tuple) and len(obj) == length def isArray(obj, shape): return isinstance(obj, numpy.ndarray) and obj.shape == shape # set up some non-trivial dimensional units agama.setUnits(length=2, velocity=3, mass=4e6) dens = agama.Density(type='plummer') pots = agama.Potential(type='dehnen', gamma=0) # spherical potential potf = agama.Potential(type='plummer', q=0.75) # flattened potential pott = agama.Potential(type='plummer', p=0.75, q=0.5) # triaxial potential actf = agama.ActionFinder(potf) actm = agama.ActionMapper(potf, [1, 1, 1]) df0 = agama.DistributionFunction(type='quasispherical', density=pots, potential=pots, r_a=2.0) df1 = agama.DistributionFunction(type='quasispherical', density=dens, potential=pots, beta0=-0.2) df2 = agama.DistributionFunction(df1, df0) # composite DF with two components gms1 = agama.GalaxyModel(pots, df1) # simple DF, spherical potential
iniDFThinDisc = dict(ini.items("DF thin disc")) iniDFThickDisc = dict(ini.items("DF thick disc")) iniDFStellarHalo = dict(ini.items("DF stellar halo")) iniDFDarkHalo = dict(ini.items("DF dark halo")) iniSCMHalo = dict(ini.items("SelfConsistentModel halo")) iniSCMDisc = dict(ini.items("SelfConsistentModel disc")) iniSCM = dict(ini.items("SelfConsistentModel")) # define external unit system describing the data (including the parameters in INI file) agama.setUnits(length=1, velocity=1, mass=1) # in Kpc, km/s, Msun # initialize the SelfConsistentModel object (only the potential expansion parameters) model = agama.SelfConsistentModel(**iniSCM) # create initial ('guessed') density profiles of all components densityBulge = agama.Density(**iniPotenBulge) densityDarkHalo = agama.Density(**iniPotenDarkHalo) densityThinDisc = agama.Density(**iniPotenThinDisc) densityThickDisc = agama.Density(**iniPotenThickDisc) densityGasDisc = agama.Density(**iniPotenGasDisc) densityStellarDisc = agama.Density(densityThinDisc, densityThickDisc) # composite # add components to SCM - at first, all of them are static density profiles model.components.append( agama.Component(dens=densityDarkHalo, disklike=False)) model.components.append( agama.Component(dens=densityStellarDisc, disklike=True)) model.components.append(agama.Component(dens=densityBulge, disklike=False)) model.components.append(agama.Component(dens=densityGasDisc, disklike=True))
) ) ### finally, decide what to do if command == 'RUN': # create a dark halo according to the provided parameters (type, scale radius and circular velocity) if rhalo>0 and vhalo>0: if halotype.upper() == 'LOG': densityHalo = agama.schwarzlib.makeDensityLogHalo(rhalo, vhalo) elif halotype.upper() == 'NFW': densityHalo = agama.schwarzlib.makeDensityNFWHalo(rhalo, vhalo) else: raise ValueError('Invalid halo type') else: densityHalo = agama.Density(type='Plummer', mass=0, scaleRadius=1) # no halo # additional density component for constructing the initial conditions: # create more orbits at small radii to better resolve the kinematics around the central black hole densityExtra = agama.Density(type='Dehnen', scaleradius=1) # fiducialMbh: Mbh used to construct initial conditions (may differ from Mbh used to integrate orbits; # the idea is to keep fiducialMbh fixed between runs with different Mbh, so that the initial conditions # for the orbit library are the same, compensating one source of noise in chi2 due to randomness) fiducialMbh = densityStars.totalMass() * 0.01 # potential of the galaxy, excluding the central BH pot_gal = agama.Potential(type='Multipole', density=agama.Density(densityStars, densityHalo), # all density components together lmax=32, mmax=0, gridSizeR=40) # mmax=0 means axisymmetry; lmax is set to a large value to accurately represent the disk # potential of the central BH
def __init__(self, filename): ''' Initialize starting values of scaled parameters and and define upper/lower limits on them; also obtain the true values by analyzing the input file name ''' self.initValues = numpy.array( [9.0, 0.1, 0.5, 4.0, 1.0, 6.0, 1.5, 1.0, 1.0, 1.0, 1.]) self.minValues = numpy.array( [6.0, -1.0, -0.5, 2.2, 0.5, 3.2, 0.0, 0.5, 0.1, 0.1, -1.]) self.maxValues = numpy.array( [10., 2.0, 1.9, 6.0, 2.5, 12., 2.8, 2.0, 2.8, 2.8, 3.]) self.labels = ( \ r'$\log_{10}(\rho_\mathrm{scale})$', r'$\log_{10}(R_\mathrm{scale})$', r'$\gamma$', r'$\beta$', r'$\alpha$', \ r'$B_\mathrm{DF}$', r'$\Gamma_\mathrm{DF}$', r'$\eta$', r'$g_r$', r'$h_r$', r'$\log_{10}(J_0)$') self.numPotParams = 5 # potential params come first in the list self.numDFParams = 6 # DF params come last in the list # parse true values m = re.search( r'gs(\d+)_bs(\d+)_rcrs([a-z\d]+)_rarc([a-z\d]+)_([a-z]+)_(\d+)mpc3', filename) n = re.search(r'data_([ch])_rh(\d+)_rs(\d+)_gs(\d+)', filename) if m: self.trueParams = ( numpy.log10(float(m.group(6)) * 1e6), 0.0, 1.0 if m.group(5) == 'cusp' else 0.0, 3.0, 1.0, # true params of tracer DF do not belong to this family of models, so ignore them numpy.nan, numpy.nan, numpy.nan, numpy.nan, numpy.nan, numpy.nan) tracerParams = dict(type='Spheroid', densityNorm=1.0, scaleRadius=float(m.group(3)) * 0.01, gamma=float(m.group(1)) * 0.01, beta=float(m.group(2)) * 0.1, alpha=2.0) beta0 = 0 r_a = float(m.group(4)) * 0.01 * float( m.group(3)) * 0.01 if m.group(4) != 'inf' else numpy.inf elif n: self.trueParams = ( numpy.log10(3.021516e7 if n.group(1) == 'c' else 2.387329e7), numpy.log10(float(n.group(2))), 0.0 if n.group(1) == 'c' else 1.0, 4.0, 1.0, numpy.nan, numpy.nan, numpy.nan, numpy.nan, numpy.nan, numpy.nan) tracerParams = dict( type='Spheroid', densityNorm=1.0, scaleRadius=1.75 if n.group(3) == '175' else 0.5, gamma=float(n.group(4)) * 0.1, beta=5.0, alpha=2.0) beta0 = -0.5 r_a = numpy.inf else: print("Can't determine true parameters!") return self.truePotential = agama.Potential( type='Spheroid', densityNorm=10**self.trueParams[0], scaleRadius=10**self.trueParams[1], gamma=self.trueParams[2], beta=self.trueParams[3], alpha=self.trueParams[4]) # normalize the tracer density profile to have a total mass of unity tracerParams["densityNorm"] /= agama.Density( **tracerParams).totalMass() self.tracerDensity = agama.Density(**tracerParams) self.tracerBeta = lambda r: (beta0 + (r / r_a)**2) / (1 + (r / r_a)**2)
def makeCMC(mass, scaleRadius, scaleHeight, axisRatioY): norm = mass / (4 * numpy.pi * scaleRadius**2 * scaleHeight * axisRatioY) return agama.Density(lambda xyz: norm * numpy.exp(-(xyz[:, 0]**2 + ( xyz[:, 1] / axisRatioY)**2)**0.5 / scaleRadius - abs(xyz[:, 2]) / scaleHeight))
# initialize the SelfConsistentModel object (only the potential expansion parameters) model = agama.SelfConsistentModel(RminCyl=0.005, RmaxCyl=1.0, sizeRadialCyl=25, zminCyl=0.005, zmaxCyl=1.0, sizeVerticalCyl=25) # construct a two component model: NSD + NSC # NSD -> generated self-consistently # NSC -> kept fixed as an external potential # NSC best-fitting model from Chatzopoulos et al. 2015 (see Equation 28 here: https://arxiv.org/pdf/2007.06577.pdf) density_NSC_init = agama.Density(type='Dehnen', mass=6.1e-3, gamma=0.71, scaleRadius=5.9e-3, axisRatioZ=0.73) # NSD model 3 from Sormani et al. 2020 (see Equation 24 here: https://arxiv.org/pdf/2007.06577.pdf) d1 = agama.Density(type='Spheroid', DensityNorm=0.9 * 222.885, gamma=0, beta=0, axisRatioZ=0.37, outerCutoffRadius=0.0050617, cutoffStrength=0.7194) d2 = agama.Density(type='Spheroid', DensityNorm=0.9 * 169.975, gamma=0, beta=0,
def __init__(self, filename): ''' Initialize starting values of scaled parameters and and define upper/lower limits on them; also obtain the true values by analyzing the input file name ''' self.initValues = numpy.array( [9.0, 0.1, 0.5, 4.0, 1.0, 0.1, 0.5, 4.0, 1.0, +0.1, 2.0]) self.minValues = numpy.array( [6.0, -1.0, -0.5, 2.2, 0.5, -2.0, 0.0, 3.5, 0.5, -0.5, -1.5]) self.maxValues = numpy.array( [10., 2.0, 1.9, 6.0, 2.5, 1.0, 1.9, 6.0, 3.0, 0.5, 5.0]) self.labels = ( \ r'$\log_{10}(\rho_\mathrm{scale})$', r'$\log_{10}(R_\mathrm{scale})$', r'$\gamma$', r'$\beta$', r'$\alpha$', \ r'$\log_{10}(R_\star)$', r'$\Gamma_\star$', r'$\mathrm{B}_\star$', r'$\mathrm{A}_\star$', r'$\beta_0$', r'$\log_{10}(r_a)$') self.numPotParams = 5 # potential params come first in the list self.numDFParams = 6 # DF params come last in the list # parse true values m = re.search( r'gs(\d+)_bs(\d+)_rcrs([a-z\d]+)_rarc([a-z\d]+)_([a-z]+)_(\d+)mpc3', filename) n = re.search(r'data_([ch])_rh(\d+)_rs(\d+)_gs(\d+)', filename) if m: self.trueParams = (numpy.log10(float(m.group(6)) * 1e6), 0.0, 1.0 if m.group(5) == 'cusp' else 0.0, 3.0, 1.0, numpy.log10(float(m.group(3)) * 0.01), float(m.group(1)) * 0.01, float(m.group(2)) * 0.10, 2.0, 0.0, numpy.log10( float(m.group(4)) * 0.01 * float(m.group(3)) * 0.01) if m.group(4) != 'inf' else 5.0) elif n: self.trueParams = ( numpy.log10(3.021516e7 if n.group(1) == 'c' else 2.387329e7), numpy.log10(float(n.group(2))), 0.0 if n.group(1) == 'c' else 1.0, 4.0, 1.0, numpy.log10(1.75 if n.group(3) == '175' else 0.5), float(n.group(4)) * 0.1, 5.0, 2.0, -0.5, 5.0) else: print("Can't determine true parameters!") return self.truePotential = agama.Potential( type='Spheroid', densityNorm=10**self.trueParams[0], scaleRadius=10**self.trueParams[1], gamma=self.trueParams[2], beta=self.trueParams[3], alpha=self.trueParams[4]) tracerParams = dict(type='Spheroid', densityNorm=1.0, scaleRadius=10**self.trueParams[self.numPotParams + 0], gamma=self.trueParams[self.numPotParams + 1], beta=self.trueParams[self.numPotParams + 2], alpha=self.trueParams[self.numPotParams + 3]) # normalize the tracer density profile to have a total mass of unity tracerParams["densityNorm"] /= agama.Density( **tracerParams).totalMass() self.tracerDensity = agama.Density(**tracerParams) self.tracerBeta = lambda r: (self.trueParams[-2] + (r / 10**self.trueParams[-1])**2) / (1 + ( r / 10**self.trueParams[-1])**2)
determined by its distribution function in terms of actions. We use the true DF of the Plummer model, expressed in terms of actions, and start iterations from a deliberately wrong initial guess; nevertheless, after 10 iterations we converge to the true solution within 1%; each iteration approximately halves the error. """ import agama, numpy, matplotlib.pyplot as plt # the distribution function defining the model truepot = agama.Potential(type='Plummer') df = agama.DistributionFunction(type='QuasiIsotropic', potential=truepot, density=truepot) # initial guess for the density profile - deliberately a wrong one dens = agama.Density(type='Dehnen', mass=0.1, scaleRadius=0.5) # define the self-consistent model consisting of a single component params = dict(rminSph=0.001, rmaxSph=1000., sizeRadialSph=40, lmaxAngularSph=0) comp = agama.Component(df=df, density=dens, disklike=False, **params) scm = agama.SelfConsistentModel(**params) scm.components = [comp] # prepare visualization r = numpy.logspace(-2., 2.) xyz = numpy.vstack((r, r * 0, r * 0)).T plt.plot(r, dens.density(xyz), label='Init density') plt.plot(r, truepot.density(xyz), label='True density', c='k')[0].set_dashes([4, 4]) # perform several iterations of self-consistent modelling procedure
def createModel(iniFileName): """ parse the ini file and construct the initial data for the entire model: all density, potentials, and Schwarzschild model components, with initial conditions for the orbit library and relevant density/kinematic targets. """ model = type('Model', (), {}) ini = RawConfigParser() ini.read(iniFileName) sec = ini.sections() sec_den = dict() # list of density components sec_pot = dict() # same for potential sec_comp = dict( ) # list of model components (each one may include several density objects) Omega = None for s in sec: if s.lower().startswith('density'): sec_den[s.lower()] = dict(ini.items(s)) if s.lower().startswith('potential'): sec_pot[s.lower()] = dict(ini.items(s)) if s.lower().startswith('component'): sec_comp[s.lower()] = dict(ini.items(s)) if s.lower() == 'global': # pattern speed for name, value in ini.items(s): if name.lower() == 'omega': Omega = float(value) # construct all density and potential objects den = dict() pot = dict() for name, value in sec_den.items(): den[name] = agama.Density(**value) for name, value in sec_pot.items(): if 'density' in value: listdens = list() for den_str in filter(None, re.split('[,\s]', value['density'].lower())): if den_str in sec_den: listdens.append(den[den_str]) if len(listdens) == 1: value['density'] = listdens[0] elif len(listdens) > 1: value['density'] = agama.Density(*listdens) # else there are no user-defined density sections print("Creating " + name) pot[name] = agama.Potential(**value) pot[name].export("model_" + name) model.density = den if len(pot) == 0: raise ValueError("No potential components defined") if len(pot) == 1: model.potential = pot.values()[0] else: model.potential = agama.Potential(*pot.values()) # determine corotation radius in case of nonzero pattern speed if Omega != 0: try: from scipy.optimize import brentq print("Omega=%.3g => corotation radius: %.3g" % ( Omega, brentq( lambda r: model.potential.force(r, 0, 0)[0] / r + Omega**2, 1e-10, 1e10, rtol=1e-4))) except Exception as e: print("Omega=%.3g, corotation radius unknown: %s" % (Omega, str(e))) # construct all model components if len(sec_comp) == 0: raise ValueError("No model components defined") model.components = dict() for name, value in sec_comp.items(): targets = list() if not 'density' in value: raise ValueError("Component " + name + " does not have an associated density") listdens = list() for den_str in filter(None, re.split('[,\s]', value['density'].lower())): if den_str in sec_den: listdens.append(den[den_str]) elif den_str in sec_pot: listdens.append(pot[den_str]) else: raise ValueError("Unknown density component: " + den_str + " in " + name) if len(listdens) == 1: density = listdens[0] elif len(listdens) > 1: density = agama.Density(*listdens) else: raise ValueError("No density components in " + name) print("Creating density target for " + name) targets.append(agama.Target(**value)) if 'kinemgrid' in value: options = { "type": 'KinemShell', \ "gridr": value['kinemgrid'], \ "degree": int(value['kinemdegree']) } print("Creating kinematic target for " + name) targets.append(agama.Target(**options)) if 'numorbits' in value: icoptions = { 'n': int(value['numorbits']), 'potential': model.potential } if 'icbeta' in value: icoptions['beta'] = float(value['icbeta']) if 'ickappa' in value: icoptions['kappa'] = float(value['ickappa']) print("Creating initial conditions for %i orbits in %s" % (icoptions['n'], name)) ic, weightprior = density.sample(**icoptions) else: raise ValueError("No orbits defined in " + name) if 'inttime' in value: inttime = float(value['inttime']) * model.potential.Tcirc(ic) else: raise ValueError("No integration time defined in " + name) comp = type('Component', (), \ {"density": density, "ic": ic, "weightprior": weightprior, \ "inttime": inttime, "targets": targets, "Omega": Omega} ) if 'trajsize' in value: comp.trajsize = int(value['trajsize']) if 'beta' in value: comp.beta = float(value['beta']) if 'nbody' in value: comp.nbody = int(value['nbody']) if not 'trajsize' in value: raise ValueError("No trajectory will be stored in " + name + ", cannot create Nbody model") model.components[name] = comp return model
def contraction(pot_dm, pot_bar, method='C20', beta_dm=0.0, rmin=1e-2, rmax=1e4): ''' Construct the contracted halo density profile for the given initial halo density and the baryonic density profiles. Arguments: pot_dm: initial halo potential (assumed to be spherical!). pot_bar: potential of baryons (a spherical approximation to it will be used even if it was not spherical). method: choice between two alternative approaches: 'C20' (default) uses the approximate correction procedure from Cautun+2020; 'adiabatic' uses the invariance of actions in conjunction with an internally constructed halo DF. beta_dm: anisotropy coefficient for the halo DF (only used with the adiabatic method, must be between -0.5 and 1, default 0 means isotropy). rmin (1e-2), rmax(1e4): min/max grid radii for representing the contracted density profile (default values are suitable for Milky Way-sized galaxies if expressed in kpc). Return: the spherically symmetric contracted halo potential. ''' gridr = numpy.logspace(numpy.log10(rmin), numpy.log10(rmax), 101) xyz = numpy.column_stack((gridr, gridr * 0, gridr * 0)) if method == 'adiabatic': # create a spherical DF for the DM-only potential/density pair with a constant anisotropy coefficient beta df_dm = agama.DistributionFunction(type='QuasiSpherical', potential=pot_dm, beta0=beta_dm) # create a sphericalized total potential (DM+baryons) pot_total_sph = agama.Potential(type='multipole', potential=agama.Potential( pot_dm, pot_bar), lmax=0, rmin=0.1 * rmin, rmax=10 * rmax) # compute the density generated by the DF in the new total potential at the radial grid dens_contracted = agama.GalaxyModel(pot_total_sph, df_dm).moments(xyz, dens=True, vel=False, vel2=False) elif method == 'C20': # use the differential (d/dr) form of Eq. (11) from Cautun et al (2020) to approximate the effect of contraction cumul_mass_dm = pot_dm.enclosedMass( gridr) # cumulative mass profile of DM cumul_mass_bar = pot_bar.enclosedMass( gridr) # same for baryons (sphericalized) valid_r = numpy.hstack([ True, cumul_mass_bar[:-1] < cumul_mass_bar[1:] * 0.999 ]) # use only those radii where mass keeps increasing sph_dens_bar = agama.Density(cumulmass=numpy.column_stack( (gridr[valid_r], cumul_mass_bar[valid_r]))) # sphericalized baryon density profile f_bar = 0.157 # cosmic baryon fraction; the formula is calibrated against simulations only for this value eta_bar = cumul_mass_bar / cumul_mass_dm * ( 1. - f_bar ) / f_bar # the last two terms account for transforming the DM mass into the corresponding baryonic mass in dark-matter-only simulations first_factor = 0.45 + 0.41 * (eta_bar + 0.98)**0.53 dens_dm_orig = pot_dm.density(xyz) temp = sph_dens_bar.density( xyz) - eta_bar * dens_dm_orig * f_bar / (1. - f_bar) const_term = 0.41 * 0.53 * (eta_bar + 0.98)**(0.53 - 1.) * ( 1. - f_bar) / f_bar * temp dens_contracted = dens_dm_orig * first_factor + const_term # new values of DM density at the radial grid else: raise RuntimeError('Invalid choice of method') # create a cubic spline interpolator in log-log space valid_r = dens_contracted > 0 # make sure the input for log-spline is positive dens_contracted_interp = agama.CubicSpline(numpy.log(gridr[valid_r]), numpy.log( dens_contracted[valid_r]), reg=True) # convert the grid-based density profile into a full-fledged potential contracted_pot = agama.Potential( type="Multipole", symmetry="spherical", rmin=rmin, rmax=rmax, density=lambda xyz: numpy.exp( dens_contracted_interp(numpy.log(numpy.sum(xyz**2, axis=1)) * 0.5) )) return contracted_pot
iniDFDarkHalo = dict(ini.items("DF dark halo")) iniDFBulge = dict(ini.items("DF bulge")) iniSCMHalo = dict(ini.items("SelfConsistentModel halo")) iniSCMBulge = dict(ini.items("SelfConsistentModel bulge")) iniSCMDisk = dict(ini.items("SelfConsistentModel disk")) iniSCM = dict(ini.items("SelfConsistentModel")) solarRadius = ini.getfloat("Data", "SolarRadius") # define external unit system describing the data (including the parameters in INI file) agama.setUnits(length=1, velocity=1, mass=1) # in Kpc, km/s, Msun # initialize the SelfConsistentModel object (only the potential expansion parameters) model = agama.SelfConsistentModel(**iniSCM) # create initial ('guessed') density profiles of all components densityBulge = agama.Density(**iniPotenBulge) densityDarkHalo = agama.Density(**iniPotenDarkHalo) densityThinDisk = agama.Density(**iniPotenThinDisk) densityThickDisk = agama.Density(**iniPotenThickDisk) densityGasDisk = agama.Density(**iniPotenGasDisk) densityStellarDisk = agama.Density(densityThinDisk, densityThickDisk) # composite # add components to SCM - at first, all of them are static density profiles model.components.append(agama.Component(density=densityStellarDisk, disklike=True)) model.components.append(agama.Component(density=densityBulge, disklike=False)) model.components.append(agama.Component(density=densityDarkHalo, disklike=False)) model.components.append(agama.Component(density=densityGasDisk, disklike=True)) # compute the initial potential model.iterate() printoutInfo(model, "init")
def createModel(iniFileName): """ parse the ini file and construct the initial data for the entire model: all density, potentials, and Schwarzschild model components, with initial conditions for the orbit library and relevant density/kinematic targets. """ ini = ConfigParser.RawConfigParser() ini.read(iniFileName) sec = ini.sections() sec_den = dict() # list of density components sec_pot = dict() # same for potential sec_comp = dict( ) # list of model components (each one may include several density objects) for s in sec: if s.lower().startswith('density'): sec_den[s.lower()] = dict(ini.items(s)) if s.lower().startswith('potential'): sec_pot[s.lower()] = dict(ini.items(s)) if s.lower().startswith('component'): sec_comp[s.lower()] = dict(ini.items(s)) # construct all density and potential objects den = dict() pot = dict() for name, value in sec_den.items(): den[name] = agama.Density(**value) for name, value in sec_pot.items(): if value.has_key('density'): listdens = list() for den_str in filter(None, re.split('[,\s]', value['density'].lower())): if sec_den.has_key(den_str): listdens.append(den[den_str]) if len(listdens) == 1: value['density'] = listdens[0] elif len(listdens) > 1: value['density'] = agama.Density(*listdens) # else there are no user-defined density sections print "Creating", name pot[name] = agama.Potential(**value) pot[name].export("model_" + name) model = {'density': den} if len(pot) == 0: raise ValueError("No potential components defined") if len(pot) == 1: model['potential'] = pot.values()[0] else: model['potential'] = agama.Potential(*pot.values()) # construct all model components if len(sec_comp) == 0: raise ValueError("No model components defined") comp = dict() model['components'] = comp for name, value in sec_comp.items(): targets = list() if not value.has_key('density'): raise ValueError("Component " + name + " does not have an associated density") listdens = list() for den_str in filter(None, re.split('[,\s]', value['density'].lower())): if sec_den.has_key(den_str): listdens.append(den[den_str]) elif sec_pot.has_key(den_str): listdens.append(pot[den_str]) else: raise ValueError("Unknown density component: " + den_str + " in " + name) if len(listdens) == 1: value['density'] = listdens[0] elif len(listdens) > 1: value['density'] = agama.Density(*listdens) else: raise ValueError("No density components in " + name) print "Creating density target for", name targets.append(agama.Target(**value)) if value.has_key('kinemshells'): options = { "type": 'KinemJeans', "gridSizeR": int(value['kinemshells']), \ "density": value['density'], "potential": model['potential'] } if value.has_key('kinemdegree'): options['degree'] = int(value['kinemdegree']) print "Creating kinematic target for", name targets.append(agama.Target(**options)) if value.has_key('numorbits'): icoptions = { 'n': int(value['numorbits']), 'potential': model['potential'] } if value.has_key('icbeta'): icoptions['beta'] = float(value['icbeta']) if value.has_key('ickappa'): icoptions['kappa'] = float(value['ickappa']) print "Creating initial conditions for", icoptions[ 'n'], "orbits in", name ic, weightprior = value['density'].sample(**icoptions) else: raise ValueError("No orbits defined in " + name) if value.has_key('inttime'): inttime = float(value['inttime']) * model['potential'].Tcirc(ic) else: raise ValueError("No integration time defined in " + name) comp[name] = { "ic": ic, "weightprior": weightprior, "inttime": inttime, "targets": targets } if value.has_key('trajsize'): comp[name]['trajsize'] = int(value['trajsize']) if value.has_key('beta'): comp[name]['beta'] = float(value['beta']) if value.has_key('nbody'): comp[name]['nbody'] = int(value['nbody']) if not value.has_key('trajsize'): raise ValueError("No trajectory will be stored in " + name + ", cannot create Nbody model") return model
ini.read(iniFileName) iniPotenHalo = dict(ini.items("Potential halo")) iniPotenBulge = dict(ini.items("Potential bulge")) iniPotenDisk = dict(ini.items("Potential disk")) iniPotenBH = dict(ini.items("Potential BH")) iniDFDisk = dict(ini.items("DF disk")) iniSCMHalo = dict(ini.items("SelfConsistentModel halo")) iniSCMBulge = dict(ini.items("SelfConsistentModel bulge")) iniSCMDisk = dict(ini.items("SelfConsistentModel disk")) iniSCM = dict(ini.items("SelfConsistentModel")) # initialize the SelfConsistentModel object (only the potential expansion parameters) model = agama.SelfConsistentModel(**iniSCM) # create initial density profiles of all components densityDisk = agama.Density(**iniPotenDisk) densityBulge = agama.Density(**iniPotenBulge) densityHalo = agama.Density(**iniPotenHalo) potentialBH = agama.Potential(**iniPotenBH) # add components to SCM - at first, all of them are static density profiles model.components.append(agama.Component(density=densityDisk, disklike=True)) model.components.append(agama.Component(density=densityBulge, disklike=False)) model.components.append(agama.Component(density=densityHalo, disklike=False)) model.components.append(agama.Component(potential=potentialBH)) # compute the initial potential model.iterate() printoutInfo(model,'init') # construct the DF of the disk component, using the initial (non-spherical) potential
#!/usr/bin/python """ Create a simple single-component spherical self-consistent model determined by its distribution function in terms of actions """ import agama, numpy, matplotlib.pyplot as plt # the distribution function defining the model df = agama.DistributionFunction(type="DoublePowerLaw", J0=1, slopeIn=0, slopeOut=6, norm=30) mass = df.totalMass() scaleradius = 2. # educated guess # initial guess for the density profile dens = agama.Density(type='Plummer', mass=mass, scaleRadius=scaleradius) print "DF mass=", mass, ' Density mass=', dens.totalMass() # should be identical # define the self-consistent model consisting of a single component params = dict(rminSph=0.01, rmaxSph=100., sizeRadialSph=25, lmaxAngularSph=0) comp = agama.Component(df=df, density=dens, disklike=False, **params) scm = agama.SelfConsistentModel(**params) scm.components=[comp] # prepare visualization r=numpy.logspace(-2.,2.) xyz=numpy.vstack((r,r*0,r*0)).T plt.plot(r, dens.density(xyz), label='Init density') # perform several iterations of self-consistent modelling procedure for i in range(5): scm.iterate()