Path: blob/main/doc/source/cookbook/power_spectrum_example.py
928 views
import matplotlib.pyplot as plt1import numpy as np23import yt45"""6Make a turbulent KE power spectrum. Since we are stratified, we use7a rho**(1/3) scaling to the velocity to get something that would8look Kolmogorov (if the turbulence were fully developed).910Ultimately, we aim to compute:11121 ^ ^*13E(k) = integral - V(k) . V(k) dS1421516n ^17where V = rho U is the density-weighted velocity field, and V is the18FFT of V.1920(Note: sometimes we normalize by 1/volume to get a spectral21energy density spectrum).222324"""252627def doit(ds):28# a FFT operates on uniformly gridded data. We'll use the yt29# covering grid for this.3031max_level = ds.index.max_level3233ref = int(np.prod(ds.ref_factors[0:max_level]))3435low = ds.domain_left_edge36dims = ds.domain_dimensions * ref3738nx, ny, nz = dims3940nindex_rho = 1.0 / 3.04142Kk = np.zeros((nx // 2 + 1, ny // 2 + 1, nz // 2 + 1))4344for vel in [("gas", "velocity_x"), ("gas", "velocity_y"), ("gas", "velocity_z")]:45Kk += 0.5 * fft_comp(46ds, ("gas", "density"), vel, nindex_rho, max_level, low, dims47)4849# wavenumbers50L = (ds.domain_right_edge - ds.domain_left_edge).d5152kx = np.fft.rfftfreq(nx) * nx / L[0]53ky = np.fft.rfftfreq(ny) * ny / L[1]54kz = np.fft.rfftfreq(nz) * nz / L[2]5556# physical limits to the wavenumbers57kmin = np.min(1.0 / L)58kmax = np.min(0.5 * dims / L)5960kbins = np.arange(kmin, kmax, kmin)61N = len(kbins)6263# bin the Fourier KE into radial kbins64kx3d, ky3d, kz3d = np.meshgrid(kx, ky, kz, indexing="ij")65k = np.sqrt(kx3d**2 + ky3d**2 + kz3d**2)6667whichbin = np.digitize(k.flat, kbins)68ncount = np.bincount(whichbin)6970E_spectrum = np.zeros(len(ncount) - 1)7172for n in range(1, len(ncount)):73E_spectrum[n - 1] = np.sum(Kk.flat[whichbin == n])7475k = 0.5 * (kbins[0 : N - 1] + kbins[1:N])76E_spectrum = E_spectrum[1:N]7778index = np.argmax(E_spectrum)79kmax = k[index]80Emax = E_spectrum[index]8182plt.loglog(k, E_spectrum)83plt.loglog(k, Emax * (k / kmax) ** (-5.0 / 3.0), ls=":", color="0.5")8485plt.xlabel(r"$k$")86plt.ylabel(r"$E(k)dk$")8788plt.savefig("spectrum.png")899091def fft_comp(ds, irho, iu, nindex_rho, level, low, delta):92cube = ds.covering_grid(level, left_edge=low, dims=delta, fields=[irho, iu])9394rho = cube[irho].d95u = cube[iu].d9697nx, ny, nz = rho.shape9899# do the FFTs -- note that since our data is real, there will be100# too much information here. fftn puts the positive freq terms in101# the first half of the axes -- that's what we keep. Our102# normalization has an '8' to account for this clipping to one103# octant.104ru = np.fft.fftn(rho**nindex_rho * u)[1050 : nx // 2 + 1, 0 : ny // 2 + 1, 0 : nz // 2 + 1106]107ru = 8.0 * ru / (nx * ny * nz)108109return np.abs(ru) ** 2110111112ds = yt.load("maestro_xrb_lores_23437")113doit(ds)114115116