Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Updates to make util_AMRgrid.py run #63

Open
wants to merge 2 commits into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
27 changes: 14 additions & 13 deletions MonteCarloMarginalizeCode/Code/RIFT/misc/amrlib.py
Original file line number Diff line number Diff line change
Expand Up @@ -301,7 +301,8 @@ def create_regular_grid_from_cell(cell, side_pts=5, return_cells=False):
# This produces a representation of the points, but not in a tuple form
grid_pts = numpy.mgrid[grid_pts]
# This tuplizes them
grid_pts = numpy.array(zip(*[grid.flatten() for grid in grid_pts]))
#grid_pts = numpy.array(zip(*[grid.flatten() for grid in grid_pts])) # python2
grid_pts = grid_pts.reshape(grid_pts.shape[0],-1).T.copy()

# useful knowledge for later
grid_spacing = numpy.array([(rb - lb) / (side_pts - 1) for lb, rb in cell._bounds])
Expand Down Expand Up @@ -377,7 +378,7 @@ def save_grid_cells_hdf(base_grp, cells, crd_sys, intr_prms=None, check=True):
assert grids.attrs["distance_coordinates"] == crd_sys

levels = []
for name, ds in grids.iteritems():
for name, ds in grids.items():
levels.append((ds.attrs["level"], name))

grid_res = numpy.diff(cells[0]._bounds).flatten()
Expand Down Expand Up @@ -429,9 +430,9 @@ def load_grid_level(h5file, level, return_labels= False, base_grp="rapidpe_grids

grids = hfile[base_grp]["grids"]
if level == -1:
level = sorted(dat.attrs["level"] for name, dat in grids.iteritems())[-1]
level = sorted(dat.attrs["level"] for name, dat in grids.items())[-1]

for name, dat in grids.iteritems():
for name, dat in grids.items():
if dat.attrs["level"] == level:
grid_res = dat.attrs["resolution"][numpy.newaxis,:]
if return_labels:
Expand Down Expand Up @@ -559,19 +560,19 @@ def check_grid(grid, intr_prms, distance_coordinates):
"""
Check the validity of points in various coordinate spaces to ensure they are physical.
"""
m1_axis, m2_axis = intr_prms.index("m1"), intr_prms.index("m2")
m1_axis, m2_axis = intr_prms.index("mass1"), intr_prms.index("mass2")
grid_check = numpy.array(grid).T
if distance_coordinates == "tau0_tau3":
bounds_mask = check_tau0tau3(grid_check[m1_axis], grid_check[m2_axis])
elif distance_coordinates == "mchirp_eta":
bounds_mask = check_mchirpeta(grid_check[m1_axis], grid_check[m2_axis])

# FIXME: Needs general spin
if "s1z" in intr_prms:
s1_axis = intr_prms.index("s1z")
if "spin1z" in intr_prms:
s1_axis = intr_prms.index("spin1z")
bounds_mask &= check_spins(grid_check[s1_axis])
if "s2z" in intr_prms:
s2_axis = intr_prms.index("s2z")
if "spin2z" in intr_prms:
s2_axis = intr_prms.index("spin2z")
bounds_mask &= check_spins(grid_check[s2_axis])
if "chi_z" in intr_prms:
chi_axis = intr_prms.index("chi_z")
Expand Down Expand Up @@ -605,10 +606,10 @@ def apply_transform(pts, intr_prms, mass_transform=None, spin_transform=None):
# You know what... recarrays are dumb, and so's your face.
# FIXME: Why does numpy want me to repack this into tuples!?
#tpts = numpy.array([tuple(pt) for pt in pts.T], dtype = numpy.dtype([(a ,"float64") for a in intr_prms]))
m1_idx, m2_idx = intr_prms.index("m1"), intr_prms.index("m2")
m1_idx, m2_idx = intr_prms.index("mass1"), intr_prms.index("mass2")
if spin_transform:
if spin_transform == "chi_z":
s1z_idx, s2z_idx = intr_prms.index("s1z"), intr_prms.index("s2z")
s1z_idx, s2z_idx = intr_prms.index("spin1z"), intr_prms.index("spin2z")
# chi_z = transform_s1zs2z_chi(pts[:,m1_idx], pts[:,m2_idx], pts[:,s1z_idx], pts[:,s2z_idx])
#The grid is converted from m1 m2 s1 s2 to mc eta chi_eff because it's better to choose the closest grid points in mc eta chi_eff space.
#chi_a is not considered when choosing the closes grid points, but we do need it to transform back to s1 s2
Expand All @@ -625,13 +626,13 @@ def apply_transform(pts, intr_prms, mass_transform=None, spin_transform=None):
return pts

def apply_inv_transform(pts, intr_prms, mass_transform=None,spin_transform=None):
m1_idx, m2_idx = intr_prms.index("m1"), intr_prms.index("m2")
m1_idx, m2_idx = intr_prms.index("mass1"), intr_prms.index("mass2")
if mass_transform:
pts[:,m1_idx], pts[:,m2_idx] = INVERSE_TRANSFORMS_MASS[VALID_TRANSFORMS_MASS[mass_transform]](pts[:,m1_idx], pts[:,m2_idx])

if spin_transform:
if spin_transform == "chi_z":
s1z_idx, s2z_idx = intr_prms.index("s1z"), intr_prms.index("s2z")
s1z_idx, s2z_idx = intr_prms.index("spin1z"), intr_prms.index("spin2z")
pts[:,s1z_idx],pts[:,s2z_idx] =transform_chi_eff_chi_a_s1zs2z(pts[:,m1_idx], pts[:,m2_idx], pts[:,s1z_idx], pts[:,s2z_idx])


Expand Down
14 changes: 7 additions & 7 deletions MonteCarloMarginalizeCode/Code/bin/util_AMRGrid.py
Original file line number Diff line number Diff line change
Expand Up @@ -155,7 +155,7 @@ def write_to_xml(cells, intr_prms, pin_prms={}, fvals=None, fname=None, verbose=
sim_insp.alpha1 = fvals[itr]
for p, v in zip(intr_prms, intr_prm._center):
setattr(sim_insp, p, v)
for p, v in pin_prms.iteritems():
for p, v in pin_prms.items():
setattr(sim_insp, p, v)
sim_insp_tbl.append(sim_insp)

Expand Down Expand Up @@ -276,12 +276,12 @@ def parse_param(popts):

#If spin 1 and 2 are not specified, they are pinned. This means the spin columns still appear in the output grid.
spin_transform=None
if not "s1z" in intr_prms or not "s2z" in intr_prms:
if not "s1z" in intr_prms and not "s2z" in intr_prms:
if not "s1z" in pin_prms:
pin_prms["s1z"] = 0.0
if not "s2z" in pin_prms:
pin_prms["s2z"] = 0.0
if not "spin1z" in intr_prms or not "spin2z" in intr_prms:
if not "spin1z" in intr_prms and not "spin2z" in intr_prms:
if not "spin1z" in pin_prms:
pin_prms["spin1z"] = 0.0
if not "spin2z" in pin_prms:
pin_prms["spin2z"] = 0.0
else:
sys.exit("spin1z or spin2z is specified but not the other spin. compute intrinsic grid is not setup to search just one")
else:
Expand Down