mirror of
https://github.com/Richard-Sti/csiborgtools.git
synced 2024-10-19 00:25:05 +02:00
e08c741fc8
* Add periodic distances * Little corrections * Fix little bug * Modernise the script * Small updates * Remove clump * Add new halo routines * Fix weights * Modernise the script * Add check ups on convergence * More convergence check ups * Edit bounds * Add default argument * Update fit heuristic and NaNs * Change maxiter * Switch NFW minimization to log-sapce * Remove print statement * Turn convert_from_box abstract property required for all boxes. * Move files * Simplify script * Improve the argument parser * Remove optinal argument * Improve argument parser * Add a minimum concentration limit
101 lines
3.5 KiB
Python
101 lines
3.5 KiB
Python
# Copyright (C) 2022 Richard Stiskalek
|
|
# This program is free software; you can redistribute it and/or modify it
|
|
# under the terms of the GNU General Public License as published by the
|
|
# Free Software Foundation; either version 3 of the License, or (at your
|
|
# option) any later version.
|
|
#
|
|
# This program is distributed in the hope that it will be useful, but
|
|
# WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General
|
|
# Public License for more details.
|
|
#
|
|
# You should have received a copy of the GNU General Public License along
|
|
# with this program; if not, write to the Free Software Foundation, Inc.,
|
|
# 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
|
|
"""
|
|
Script to sort the initial snapshot particles according to their final
|
|
snapshot ordering, which is sorted by the halo IDs.
|
|
"""
|
|
from argparse import ArgumentParser
|
|
from datetime import datetime
|
|
from gc import collect
|
|
|
|
import h5py
|
|
import numpy
|
|
from mpi4py import MPI
|
|
from taskmaster import work_delegation
|
|
|
|
from utils import get_nsims
|
|
|
|
try:
|
|
import csiborgtools
|
|
except ModuleNotFoundError:
|
|
import sys
|
|
|
|
sys.path.append("../")
|
|
import csiborgtools
|
|
|
|
|
|
def _main(nsim, simname, verbose):
|
|
"""
|
|
Sort the initial snapshot particles according to their final snapshot
|
|
ordering and dump them into a HDF5 file.
|
|
|
|
Parameters
|
|
----------
|
|
nsim : int
|
|
IC realisation index.
|
|
simname : str
|
|
Simulation name.
|
|
verbose : bool
|
|
Verbosity flag.
|
|
"""
|
|
if simname == "quijote":
|
|
raise NotImplementedError("Quijote not implemented yet.")
|
|
|
|
paths = csiborgtools.read.Paths(**csiborgtools.paths_glamdring)
|
|
partreader = csiborgtools.read.ParticleReader(paths)
|
|
|
|
if verbose:
|
|
print(f"{datetime.now()}: reading and processing simulation {nsim}.",
|
|
flush=True)
|
|
# We first load the particle IDs in the final snapshot.
|
|
pidf = csiborgtools.read.read_h5(paths.particles(nsim))
|
|
pidf = pidf["particle_ids"]
|
|
# Then we load the particles in the initil snapshot and make sure that
|
|
# their particle IDs are sorted as in the final snapshot. Again, because of
|
|
# precision this must be read as structured.
|
|
# NOTE: ID has to be the last column.
|
|
pars_extract = ["x", "y", "z", "M", "ID"]
|
|
part0, pid0 = partreader.read_particle(
|
|
1, nsim, pars_extract, return_structured=False, verbose=verbose)
|
|
# First enforce them to already be sorted and then apply reverse
|
|
# sorting from the final snapshot.
|
|
part0 = part0[numpy.argsort(pid0)]
|
|
del pid0
|
|
collect()
|
|
part0 = part0[numpy.argsort(numpy.argsort(pidf))]
|
|
if verbose:
|
|
print(f"{datetime.now()}: dumping particles for {nsim}.", flush=True)
|
|
with h5py.File(paths.initmatch(nsim, "particles"), "w") as f:
|
|
f.create_dataset("particles", data=part0)
|
|
|
|
|
|
if __name__ == "__main__":
|
|
# Argument parser
|
|
parser = ArgumentParser()
|
|
parser.add_argument("--simname", type=str, default="csiborg",
|
|
choices=["csiborg", "quijote"],
|
|
help="Simulation name")
|
|
parser.add_argument("--nsims", type=int, nargs="+", default=None,
|
|
help="IC realisations. If `-1` processes all.")
|
|
args = parser.parse_args()
|
|
|
|
paths = csiborgtools.read.Paths(**csiborgtools.paths_glamdring)
|
|
nsims = get_nsims(args, paths)
|
|
|
|
def main(nsim):
|
|
_main(nsim, args.simname, MPI.COMM_WORLD.Get_size() == 1)
|
|
|
|
work_delegation(main, nsims, MPI.COMM_WORLD)
|