Skip to content
Merged
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
32 changes: 32 additions & 0 deletions scri/asymptotic_bondi_data/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -175,6 +175,38 @@ def psi0(self, psi0prm):
self._psi0[:] = psi0prm
return self.psi0

#Slicing
def __getitem__(self, key):
"""
Extract time slices of the asymptotic Bondi data efficiently.
"""
# If key is a valid time slice or index, extract the corresponding
# sliced data

if not isinstance(key, (slice, int)):
raise ValueError(f"Invalid key `{key}` of type `{type(key)}`.")

import copy
import functools

ModesTS = functools.partial(ModesTimeSeries, ell_max=self.ell_max, multiplication_truncator=self.sigma._metadata['multiplication_truncator'])

new_abd = copy.copy(self)
new_abd._raw_data = self._raw_data[:, key, :]
new_abd._time = self._time[key]

new_abd._psi0 = ModesTS(new_abd._raw_data[0], new_abd._time, spin_weight=2)
new_abd._psi1 = ModesTS(new_abd._raw_data[1], new_abd._time, spin_weight=1)
new_abd._psi2 = ModesTS(new_abd._raw_data[2], new_abd._time, spin_weight=0)
new_abd._psi3 = ModesTS(new_abd._raw_data[3], new_abd._time, spin_weight=-1)
new_abd._psi4 = ModesTS(new_abd._raw_data[4], new_abd._time, spin_weight=-2)
new_abd._sigma = ModesTS(new_abd._raw_data[5], new_abd._time, spin_weight=2)

if self.frame.shape[0] == self.n_times:
new_abd.frame = self.frame[key]
return new_abd


def copy(self):
import copy

Expand Down
23 changes: 11 additions & 12 deletions scri/asymptotic_bondi_data/map_to_superrest_frame.py
Original file line number Diff line number Diff line change
Expand Up @@ -783,9 +783,8 @@ def map_to_superrest_frame(
if alpha_ell_max is None:
alpha_ell_max = ell_max

abd_interp = abd.interpolate(
abd.t[np.argmin(abs(abd.t - (t_0 - (padding_time + 200)))) : np.argmin(abs(abd.t - (t_0 + (padding_time + 200)))) + 1]
)
# slice the abd object using padding time
abd_sliced = abd[np.abs(abd.t - (t_0 - (padding_time + 200))).argmin() : np.abs(abd.t - (t_0 + (padding_time + 200))).argmin() + 1]

# apply a time translation so that we're mapping
# to the superrest frame at u = 0
Expand All @@ -810,7 +809,7 @@ def map_to_superrest_frame(
pass

if itr == 0:
abd_interp_prime = abd_interp.transform(
abd_sliced_prime = abd_sliced.transform(
supertranslation=BMS_transformation.supertranslation,
frame_rotation=BMS_transformation.frame_rotation.components,
boost_velocity=BMS_transformation.boost_velocity,
Expand All @@ -819,7 +818,7 @@ def map_to_superrest_frame(
for transformation in order:
if transformation == "supertranslation":
new_transformation, supertranslation_rel_errs = supertranslation_to_map_to_superrest_frame(
abd_interp_prime,
abd_sliced_prime,
target_PsiM,
N_itr_max=N_itr_maxes["supertranslation"],
rel_err_tol=rel_err_tols["supertranslation"],
Expand All @@ -828,27 +827,27 @@ def map_to_superrest_frame(
)
elif transformation == "rotation":
new_transformation, rot_rel_errs = rotation_to_map_to_superrest_frame(
abd_interp_prime,
abd_sliced_prime,
target_strain=target_strain,
N_itr_max=N_itr_maxes["rotation"],
rel_err_tol=rel_err_tols["rotation"],
print_conv=print_conv,
)
elif transformation == "CoM_transformation":
new_transformation, CoM_rel_errs = com_transformation_to_map_to_superrest_frame(
abd_interp_prime,
abd_sliced_prime,
N_itr_max=N_itr_maxes["CoM_transformation"],
rel_err_tol=rel_err_tols["CoM_transformation"],
print_conv=print_conv,
)
elif transformation == "time_phase":
if target_strain is not None:
strain_interp_prime = scri.asymptotic_bondi_data.map_to_superrest_frame.MT_to_WM(
2.0 * abd_interp_prime.sigma.bar, dataType=scri.h
strain_sliced_prime = scri.asymptotic_bondi_data.map_to_superrest_frame.MT_to_WM(
2.0 * abd_sliced_prime.sigma.bar, dataType=scri.h
)

rel_err, _, res = align2d(
MT_to_WM(WM_to_MT(strain_interp_prime), sxs_version=True),
MT_to_WM(WM_to_MT(strain_sliced_prime), sxs_version=True),
MT_to_WM(WM_to_MT(target_strain), sxs_version=True),
0 - padding_time,
0 + padding_time,
Expand All @@ -867,7 +866,7 @@ def map_to_superrest_frame(
["supertranslation", "frame_rotation", "boost_velocity"]
)

abd_interp_prime = abd_interp.transform(
abd_sliced_prime = abd_sliced.transform(
supertranslation=BMS_transformation.supertranslation,
frame_rotation=BMS_transformation.frame_rotation.components,
boost_velocity=BMS_transformation.boost_velocity,
Expand All @@ -877,7 +876,7 @@ def map_to_superrest_frame(
# rel_err is obtained from align2d, so do nothing
pass
else:
rel_err = rel_err_for_abd_in_superrest(abd_interp_prime, target_PsiM, target_strain)
rel_err = rel_err_for_abd_in_superrest(abd_sliced_prime, target_PsiM, target_strain)

if np.mean(rel_err) < min([np.mean(r) for r in rel_errs]):
best_BMS_transformation = BMS_transformation.copy()
Expand Down