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

add show_wfs_during_program function #798

Merged
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
95 changes: 95 additions & 0 deletions webbpsf/mast_wss.py
Original file line number Diff line number Diff line change
Expand Up @@ -221,6 +221,7 @@
return mast_retrieve_opd(closest_fn, output_path = output_path)



### Functions for format conversion of OPDs

def import_wss_opd(filename, npix_out=1024, verbose=False):
Expand Down Expand Up @@ -621,3 +622,97 @@
ta_hdul = fits.open(mast_file_url)

return ta_hdul


# Functions for retrieving metadata about science observations

def _query_program_visit_times_by_inst(program, instrument, verbose=False):
""" Get the start and end times of all completed visits in a program, per instrument.
Not intended for general use; this is mostly a helper to query_program_visit_times.

Getting the vststart_mjd and visitend_mjd fields requires using the instrument keywords
interface, so one has to specify which instrument ahead of time.

Parameters
----------
program : int or str
Program ID
instrument : str
instrument name
verbose : bool
be more verbose in output?

returns list of (visitid, start, end) tuples.

"""

from astroquery.mast import Mast
svc_table = {'MIRI':'Mast.Jwst.Filtered.Miri',

Check warning on line 650 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L649-L650

Added lines #L649 - L650 were not covered by tests
'NIRCAM': 'Mast.Jwst.Filtered.NIRCam',
'NIRSPEC': 'Mast.Jwst.Filtered.NIRSpec',
'NIRISS': 'Mast.Jwst.Filtered.NIRISS',
}

service = svc_table[instrument.upper()]

Check warning on line 656 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L656

Added line #L656 was not covered by tests

collist = 'filename, program, observtn, visit_id, vststart_mjd, visitend_mjd, bstrtime'
all_columns = False

Check warning on line 659 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L658-L659

Added lines #L658 - L659 were not covered by tests

def set_params(parameters):
return [{"paramName" : p, "values" : v} for p, v in parameters.items()]

Check warning on line 662 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L661-L662

Added lines #L661 - L662 were not covered by tests


keywords = {'program': [str(program),]}
parameters = {'columns': '*' if all_columns else collist,

Check warning on line 666 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L665-L666

Added lines #L665 - L666 were not covered by tests
'filters': set_params(keywords)}

if verbose:
print("MAST query parameters:")
print(parameters)

Check warning on line 671 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L669-L671

Added lines #L669 - L671 were not covered by tests

responsetable = Mast.service_request(service, parameters)
responsetable.sort(keys='bstrtime')

Check warning on line 674 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L673-L674

Added lines #L673 - L674 were not covered by tests


visit_times = []

Check warning on line 677 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L677

Added line #L677 was not covered by tests

for row in responsetable:
visit_times.append( ('V'+row['visit_id'], row['vststart_mjd'], row['visitend_mjd']))

Check warning on line 680 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L679-L680

Added lines #L679 - L680 were not covered by tests

visit_times= set(visit_times)
return list(visit_times)

Check warning on line 683 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L682-L683

Added lines #L682 - L683 were not covered by tests

def query_program_visit_times(program, verbose=False):
""" Get the start and end times of all completed visits in a program.

Parameters
----------
program : int or str
Program ID
verbose : bool
be more verbose in output?

Returns astropy Table with columns for visit ID and start and end times.
"""

from astroquery.mast import Observations
obs = Observations.query_criteria(obs_collection=["JWST"], proposal_id=[program])

Check warning on line 699 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L698-L699

Added lines #L698 - L699 were not covered by tests
# Annoyingly, that query interface doesn't return start/end times
instruments = set([val.split('/')[0] for val in set(obs['instrument_name'])])

Check warning on line 701 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L701

Added line #L701 was not covered by tests

visit_times = []
for inst in instruments:
if verbose:
print(f"querying for visits using {inst}")
visit_times += _query_program_visit_times_by_inst(program, inst)

Check warning on line 707 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L703-L707

Added lines #L703 - L707 were not covered by tests


vids = [v[0] for v in visit_times]
starts =astropy.time.Time([float(v[1]) for v in visit_times], format='mjd')
ends = astropy.time.Time([float(v[2]) for v in visit_times], format='mjd')

Check warning on line 712 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L710-L712

Added lines #L710 - L712 were not covered by tests

#visit_times = np.asarray(visit_times)
return astropy.table.Table([vids, starts, ends],

Check warning on line 715 in webbpsf/mast_wss.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/mast_wss.py#L715

Added line #L715 was not covered by tests
names=('visit_id', 'start_mjd', 'end_mjd'))


168 changes: 162 additions & 6 deletions webbpsf/trending.py
Original file line number Diff line number Diff line change
Expand Up @@ -952,12 +952,11 @@
return start_date, end_date


def filter_opdtable_for_month(year, mon, opdtable):
"""Filter existing opdtable for a given month
This includes the last measurement in the prior month too (if applicable), so we can compute a delta
def filter_opdtable_for_daterange(start_date, end_date, opdtable):
"""Filter existing opdtable for a given time range
This includes the last measurement in the prior time range too (if applicable), so we can compute a delta
to the first one
"""
start_date, end_date = get_month_start_end(year, mon)
# Start a little early, such that we are going to have at least 1 WFS before the start date
pre_start_date = astropy.time.Time(start_date) - astropy.time.TimeDelta(4 * u.day)
opdtable = webbpsf.mast_wss.filter_opd_table(
Expand All @@ -975,6 +974,30 @@
return opdtable


def filter_opdtable_for_month(year, mon, opdtable):
"""Filter existing opdtable for a given month
This includes the last measurement in the prior month too (if applicable), so we can compute a delta
to the first one
"""
start_date, end_date = get_month_start_end(year, mon)
return filter_opdtable_for_daterange(start_date, end_date, opdtable)



def get_opdtable_for_daterange(start_date, end_date):
"""Return table of OPD measurements for date range.

This includes the last measurement preceding this date range, too, so we
can compute the first delta at the start of this range.
"""
# Retrieve full OPD table, then trim to the selected time period
opdtable0 = webbpsf.mast_wss.retrieve_mast_opd_table()
opdtable0 = webbpsf.mast_wss.deduplicate_opd_table(opdtable0)

Check warning on line 995 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L994-L995

Added lines #L994 - L995 were not covered by tests

opdtable = filter_opdtable_for_daterange(start_date, end_date, opdtable0)
return opdtable

Check warning on line 998 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L997-L998

Added lines #L997 - L998 were not covered by tests


def get_opdtable_for_month(year, mon):
"""Return table of OPD measurements for a given month.
retrieve the opd table and filter according to month/year designation
Expand Down Expand Up @@ -1043,8 +1066,6 @@
return




def monthly_trending_plot(year, month, verbose=True, instrument='NIRCam', filter='F200W', vmax=200, pid=None, opdtable=None):
"""Make monthly trending plot showing OPDs, mirror moves, RMS WFE, and the resulting PSF EEs

Expand Down Expand Up @@ -1416,6 +1437,7 @@

return fig


def plot_wfs_obs_delta(fn1, fn2, vmax_fraction=1.0):
""" Display comparison of two weak lens observations

Expand Down Expand Up @@ -1516,6 +1538,8 @@
"""Make a helpful plot showing available WFS before and after some given science
observation. This can be used to help inform how much WFE variability there was around that time.

See also show_wfs_during_program.

Parameters
----------
filename : str
Expand Down Expand Up @@ -1624,6 +1648,138 @@
ax4.set_title("Delta WFE\nAfter-Before", color='C1', fontweight='bold')


def show_wfs_during_program(program, verbose=False, ax = None, ref_wavefront_date=None,
ref_wavefront_visit=None,
start_date=None, end_date=None):
""" Show WFS data for the entire time interval in which a given program was observed.

Plots time series of the WFS measured RMS WFE as seen in NIRCam, the start and end times
of all the observations in that program, and the delta wavefront RMS relative to either
the median wavefront during that time period, or to a specifed date.

See also show_wfs_around_obs.

Parameters
----------
program : str or int
Program ID number
verbose : str
be more verbose in output?
ax : matplotlib.Axes instance, or None
an existing Axes to plot into, or else a new one will be created.
ref_wavefront_date : date-like str or None
Optional, to specify which date's wavefront sensing should be used as
the reference wavefront for computation of the plotted delta WFE RMS.
The closest wavefront in time to the specified date will be used.
If not set, the median wavefront over the entire time period will be used.
ref_wavefront_visit : str
Another way to specify which date's wavefront sensing should be used as the
reference. Provide a science visit ID from the program (e.g. "V01234002001")
and the closest wavefront to the observation start time of that visit will be used.
start_date, end_date : strings or astropy.time.Time
Start and/or end dates for the time period to display. If not set,
reasonable default values will be computed that include all observations
for that science program plus some padding time on either side.
"""

# Query mast for when the observations took place
science_visit_table = webbpsf.mast_wss.query_program_visit_times(program, verbose=verbose)
science_visit_table.sort(keys=['start_mjd'])
if verbose:
for row in science_visit_table:
print(f" Found {row['visit_id']} starting at {row['start_mjd'].iso}")

Check warning on line 1690 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1686-L1690

Added lines #L1686 - L1690 were not covered by tests


# Figure out reasonable start and end dates for the time interval to display,
# or use values provided by the user
sci_duration = science_visit_table['start_mjd'].max() - science_visit_table['start_mjd'].min()
plot_padding_time_range = max(sci_duration.to(u.day)*0.2, 4*u.day)
if start_date is None:
start_date = science_visit_table['start_mjd'].min() - plot_padding_time_range

Check warning on line 1698 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1695-L1698

Added lines #L1695 - L1698 were not covered by tests
else:
start_date = astropy.time.Time(start_date)
if end_date is None:
end_date = science_visit_table['end_mjd'].max() + plot_padding_time_range

Check warning on line 1702 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1700-L1702

Added lines #L1700 - L1702 were not covered by tests
else:
end_date = astropy.time.Time(end_date)

Check warning on line 1704 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1704

Added line #L1704 was not covered by tests

# Look up wavefront sensing and mirror move corrections for that range
opdtable = get_opdtable_for_daterange(start_date, end_date)
corrections_table = webbpsf.mast_wss.get_corrections(opdtable)

Check warning on line 1708 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1707-L1708

Added lines #L1707 - L1708 were not covered by tests

# Iterate over the WFS measurements to retrieve the OPDs and RMS WFE
wfs_dates = []
rms_obs = []
opds = []
for row_index in range(len(opdtable)):
opd_fn = opdtable[row_index]['fileName']
try:
opd, opd_hdul = webbpsf.trending._read_opd(opd_fn)
except FileNotFoundError:
webbpsf.mast_wss.mast_retrieve_opd(opd_fn, verbose=verbose)
opd, opd_hdul = webbpsf.trending._read_opd(opd_fn)

Check warning on line 1720 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1711-L1720

Added lines #L1711 - L1720 were not covered by tests

opds.append(opd)
wfs_dates.append(opdtable[row_index]['date'])
rms_obs.append(opd_hdul[1].header['RMS_WFE']*1000)
wfs_dates_array = astropy.time.Time(wfs_dates, format='isot')
opd_array = np.asarray(opds)

Check warning on line 1726 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1722-L1726

Added lines #L1722 - L1726 were not covered by tests

# Compute delta WFE, either relative to median or a specifed date
if ref_wavefront_visit is not None:
for row in science_visit_table:
if row['visit_id']==ref_wavefront_visit:
ref_wavefront_date = row['start_mjd'].iso
if verbose:
print(f" Will compare wavefronts relative to {row['visit_id']} starting at {row['start_mjd'].iso}")

Check warning on line 1734 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1729-L1734

Added lines #L1729 - L1734 were not covered by tests


if ref_wavefront_date:
refdate = astropy.time.Time(ref_wavefront_date)
delta_times = wfs_dates_array - refdate
closest = np.argmin(np.abs(delta_times))
reference_opd = opds[closest]
ref_label = "WFS near "+ref_wavefront_date
if verbose:
print(f"Computing delta OPDs relative to date {ref_wavefront_date}, using opd {opdtable[closest]['fileName']}")

Check warning on line 1744 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1737-L1744

Added lines #L1737 - L1744 were not covered by tests
else:
median_opd = np.median(opd_array, axis=0)
reference_opd = median_opd
ref_label = 'median wavefront'
if verbose:
print("Computing delta OPDs relative to median OPD over that time period.")

Check warning on line 1750 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1746-L1750

Added lines #L1746 - L1750 were not covered by tests

delta_opds = opd_array - reference_opd
mask = opd_array.sum(axis=0)!=0
delta_rmses = [webbpsf.utils.rms(d, mask=mask)*1000 for d in delta_opds]

Check warning on line 1754 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1752-L1754

Added lines #L1752 - L1754 were not covered by tests

# Plot!
if ax is None:
fig, ax = plt.subplots(figsize=(8,4), ncols=1, nrows=1)

Check warning on line 1758 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1757-L1758

Added lines #L1757 - L1758 were not covered by tests

ax.plot_date(wfs_dates_array.plot_date, rms_obs, '+',

Check warning on line 1760 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1760

Added line #L1760 was not covered by tests
color='C1', ls='-', label='Measured RMS Wavefront Error at NIRCam NRCA3')
ax.plot_date(wfs_dates_array.plot_date, delta_rmses,'none',

Check warning on line 1762 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1762

Added line #L1762 was not covered by tests
color='C0', ls='--', label=f'RMS Delta Wavefront relative to {ref_label}')

plot_sci_y = 25
ax.scatter(science_visit_table['start_mjd'].plot_date,

Check warning on line 1766 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1765-L1766

Added lines #L1765 - L1766 were not covered by tests
[plot_sci_y]*len(science_visit_table), s=150, marker="*", color='black' ,
label=f'Program {program} observations')
for row in science_visit_table:
ax.fill_betweenx([0,120], row['start_mjd'].plot_date, (row['end_mjd']).plot_date,

Check warning on line 1770 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1769-L1770

Added lines #L1769 - L1770 were not covered by tests
color='gray', alpha=0.2)
ax.text(row['start_mjd'].plot_date, plot_sci_y+4, row['visit_id'], rotation=45, fontsize='small')

Check warning on line 1772 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1772

Added line #L1772 was not covered by tests

ax.set_ylim(0,120)
ax.legend(framealpha=0.99)
ax.xaxis.set_major_formatter(

Check warning on line 1776 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1774-L1776

Added lines #L1774 - L1776 were not covered by tests
matplotlib.dates.ConciseDateFormatter(ax.xaxis.get_major_locator()))
ax.set_ylabel("Measured Wavefront Error RMS [nm]")
ax.set_xlim(start_date.plot_date, end_date.plot_date)

Check warning on line 1779 in webbpsf/trending.py

View check run for this annotation

Codecov / codecov/patch

webbpsf/trending.py#L1778-L1779

Added lines #L1778 - L1779 were not covered by tests



#### Functions for image comparisons
def show_wfs_ta_img(visitid, ax=None, return_handles=False):
""" Retrieve and display a WFS target acq image"""
Expand Down
Loading