Add CF4 mean & std (#132)

* Add observer velocity

* Add radec to supergalactic

* Update script

* Add CF4 params

* Add paths to CF4

* Add CF4 but not sth wrong

* Update nb

* Remove errors

* Update nb

* Update scripts
This commit is contained in:
Richard Stiskalek 2024-06-27 15:02:33 +01:00 committed by GitHub
parent 76e8f71556
commit 301ff61e89
No known key found for this signature in database
GPG Key ID: B5690EEEBB952194
9 changed files with 284 additions and 110 deletions

View File

@ -41,7 +41,7 @@ from sklearn.model_selection import KFold
from tqdm import trange
from ..params import SPEED_OF_LIGHT, simname2Omega_m
from ..utils import fprint, radec_to_galactic
from ..utils import fprint, radec_to_galactic, radec_to_supergalactic
H0 = 100 # km / s / Mpc
@ -105,7 +105,11 @@ class DataLoader:
# In case of Carrick 2015 the box is in galactic coordinates..
if simname == "Carrick2015":
# Carrick+2015 box is in galactic coordinates
d1, d2 = radec_to_galactic(self._cat["RA"], self._cat["DEC"])
elif "CF4" in simname:
# CF4 box is in supergalactic coordinates
d1, d2 = radec_to_supergalactic(self._cat["RA"], self._cat["DEC"])
else:
d1, d2 = self._cat["RA"], self._cat["DEC"]
@ -129,10 +133,15 @@ class DataLoader:
mean_rho_matter *= self._Omega_m
self._los_density /= mean_rho_matter
# Since Carrick+2015 provide `rho / <rho> - 1`
if simname == "Carrick2015":
# Since Carrick+2015 and CF4 provide `rho / <rho> - 1`
if simname in ["Carrick2015", "CF4", "CF4gp"]:
self._los_density += 1
# But some CF4 delta values are < -1. Check that CF4 really reports
# this.
if simname in ["CF4", "CF4gp"]:
self._los_density = np.clip(self._los_density, 1e-5, None,)
self._mask = np.ones(len(self._cat), dtype=bool)
self._catname = catalogue

View File

@ -77,6 +77,8 @@ def simname2boxsize(simname):
"quijote": 1000.,
"TNG300-1": 205.,
"Carrick2015": 400.,
"CF4": 1000., # These need to be checked with Helene Courtois.
"CF4gp": 1000.,
}
boxsize = d.get(simname, None)
@ -97,6 +99,8 @@ def simname2Omega_m(simname):
"borg2": 0.3111,
"borg2_all": 0.3111,
"Carrick2015": 0.3,
"CF4": 0.3,
"CF4gp": 0.3,
}
omega_m = d.get(simname, None)

View File

@ -115,6 +115,8 @@ class Paths:
files = [int(search(r'chain_(\d+)', f).group(1)) for f in files]
elif simname == "Carrick2015":
return [0]
elif simname in ["CF4", "CF4gp"]:
return [0]
else:
raise ValueError(f"Unknown simulation name `{simname}`.")

View File

@ -163,6 +163,12 @@ def radec_to_galactic(ra, dec):
return c.galactic.l.degree, c.galactic.b.degree
def radec_to_supergalactic(ra, dec):
"""Convert right ascension and declination to supergalactic coordinates."""
c = SkyCoord(ra=ra*u.degree, dec=dec*u.degree, frame='icrs')
return c.supergalactic.sgl.degree, c.supergalactic.sgb.degree
@jit(nopython=True, fastmath=True, boundscheck=False)
def great_circle_distance(x1, x2, in_degrees=True):
"""

File diff suppressed because one or more lines are too long

View File

@ -28,7 +28,8 @@ from astropy import units as u
from astropy.coordinates import SkyCoord
from h5py import File
from mpi4py import MPI
from taskmaster import work_delegation
from taskmaster import work_delegation # noqa
from astropy.io import fits
from utils import get_nsims
@ -88,14 +89,17 @@ def get_los(catalogue_name, simname, comm):
if comm.Get_rank() == 0:
print(f"The dataset contains {len(RA)} objects.")
# The Carrick+2015 is in galactic coordinates, so we need to convert
# the RA/dec to galactic coordinates.
if simname == "Carrick2015":
# The Carrick+2015 is in galactic coordinates, so we need to
# convert the RA/dec to galactic coordinates.
c = SkyCoord(ra=RA*u.degree, dec=dec*u.degree, frame='icrs')
pos = np.vstack((c.galactic.l, c.galactic.b)).T
elif "CF4" in simname:
# CF4 fields are in supergalactic coordinates.
c = SkyCoord(ra=RA*u.degree, dec=dec*u.degree, frame='icrs')
pos = np.vstack((c.supergalactic.sgl, c.supergalactic.sgb)).T
else:
pos = np.vstack((RA, dec)).T
else:
pos = None
@ -155,6 +159,25 @@ def get_field(simname, nsim, kind, MAS, grid):
return field
else:
raise ValueError(f"Unknown field kind: `{kind}`.")
elif "CF4" in simname:
folder = "/mnt/extraspace/rstiskalek/catalogs"
warn(f"Using local paths from `{folder}`.", RuntimeWarning)
if kind == "density":
fpath = join(folder, "CF4_new_64-z008_delta.fits")
elif kind == "velocity":
fpath = join(folder, "CF4_new_64-z008_velocity.fits")
else:
raise ValueError(f"Unknown field kind: `{kind}`.")
fpath = fpath.replace("CF4", "CF4gp") if "CF4gp" in simname else fpath
field = fits.open(fpath)[0].data
# https://projets.ip2i.in2p3.fr//cosmicflows/ says to multiply by 52
if kind == "velocity":
field *= 52
return field.astype(np.float32)
else:
raise ValueError(f"Unknown simulation name: `{simname}`.")

View File

@ -1,6 +1,6 @@
nthreads=1
memory=7
on_login=0
on_login=1
queue="berg"
env="/mnt/users/rstiskalek/csiborgtools/venv_csiborg/bin/python"
file="field_los.py"

View File

@ -17,11 +17,11 @@ if [ "$on_login" != "1" ] && [ "$on_login" != "0" ]; then
fi
# Submit a job for each combination of simname, catalogue, ksim
for simname in "csiborg2_main"; do
for catalogue in "2MTF"; do
for simname in "CF4gp"; do
for catalogue in "LOSS"; do
# for ksim in 0 1 2; do
for ksim in 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 "none"; do
# for ksim in 0; do
# for ksim in 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 "none"; do
for ksim in 0; do
pythoncm="$env $file --catalogue $catalogue --simname $simname --ksim $ksim --ksmooth $ksmooth --ndevice $ndevice --device $device"
if [ $on_login -eq 1 ]; then

View File

@ -89,6 +89,7 @@ def main(nsim, folder, fname_basis, Rmax, subtract_observer_velocity,
bf_vrad_weighted_part = np.full_like(bf_volume_part, np.nan)
bf_vrad_weighted_halo_uniform = np.full_like(bf_volume_part, np.nan)
bf_vrad_weighted_halo = np.full_like(bf_volume_part, np.nan)
obs_vel = np.full((len(observers), 3), np.nan)
for i in range(len(observers)):
print(f"{t()}: Calculating bulk flow for observer {i + 1} of simulation {nsim}.") # noqa
@ -115,13 +116,17 @@ def main(nsim, folder, fname_basis, Rmax, subtract_observer_velocity,
halo_mass_current = halo_mass[indxs]
# Subtract the observer velocity
if subtract_observer_velocity:
rscale = 0.5 # Mpc / h
weights = np.exp(-0.5 * (rdist_part / rscale)**2)
obs_vel_x = np.average(part_vel_current[:, 0], weights=weights)
obs_vel_y = np.average(part_vel_current[:, 1], weights=weights)
obs_vel_z = np.average(part_vel_current[:, 2], weights=weights)
rscale = 0.5 # Mpc / h
weights = np.exp(-0.5 * (rdist_part / rscale)**2)
obs_vel_x = np.average(part_vel_current[:, 0], weights=weights)
obs_vel_y = np.average(part_vel_current[:, 1], weights=weights)
obs_vel_z = np.average(part_vel_current[:, 2], weights=weights)
obs_vel[i, 0] = obs_vel_x
obs_vel[i, 1] = obs_vel_y
obs_vel[i, 2] = obs_vel_z
if subtract_observer_velocity:
part_vel_current[:, 0] -= obs_vel_x
part_vel_current[:, 1] -= obs_vel_y
part_vel_current[:, 2] -= obs_vel_z
@ -168,6 +173,7 @@ def main(nsim, folder, fname_basis, Rmax, subtract_observer_velocity,
f["bf_volume_halo_uniform"] = bf_volume_halo_uniform
f["bf_vrad_weighted_halo_uniform"] = bf_vrad_weighted_halo_uniform
f["bf_vrad_weighted_halo"] = bf_vrad_weighted_halo
f["obs_vel"] = obs_vel
for i in range(len(observers)):
g = f.create_group(f"obs_{str(i)}")