#!/usr/bin/env python
#
# FOLLOWUP TEST (if known injection)
#    util_PrintInjection.py --inj mdc.xml.gz --event 0 --verbose
#   util_PrintInjection.py --inj maxpt_zero_noise.xml.gz.xml.gz --event 0
"""
Integrate the extrinsic parameters of the prefactored likelihood function.
"""

import sys
import functools
from optparse import OptionParser, OptionGroup

import numpy
import numpy as np
try:
  import cupy
  xpy_default=cupy
  identity_convert = cupy.asnumpy
  junk_to_check_installed = cupy.array(5)  # this will fail if GPU not installed correctly

  print(cupy.show_config())  # print provenance/debugging information

  # Check memory allocation
  mem_info = cupy.cuda.Device().mem_info  # memory available in bytes
  n_mb_total = mem_info[1]/1024./1024.
  n_mb_used = mem_info[0]/1024/1024.
  print(  " cupy memory [total, available] {} {} Mb".format(n_mb_total,n_mb_total-n_mb_used))
  # Add failure mode test if memory insufficient!
  if n_mb_total - n_mb_used < 3000:
    print( "  - low cupy memory - ")
  cupy_success=True
except:
  print( ' no cupy')
#  import numpy as cupy  # will automatically replace cupy calls with numpy!
  xpy_default=numpy  # just in case, to make replacement clear and to enable override
  identity_convert = lambda x: x  # trivial return itself
  cupy_success=False

import lal
from glue.ligolw import utils, lsctables, table, ligolw
from glue.ligolw.utils import process
import glue.lal

import RIFT.lalsimutils as lalsimutils
import RIFT.likelihood.factored_likelihood as factored_likelihood
import RIFT.integrators.mcsampler as mcsampler
try:
    import RIFT.integrators.mcsamplerEnsemble as mcsamplerEnsemble
    mcsampler_gmm_ok = True
except:
    print(" No mcsamplerEnsemble ")
    mcsampler_gmm_ok = False
try:
    import RIFT.integrators.mcsamplerGPU as mcsamplerGPU
    mcsampler_gpu_ok = True
except:
    print( " No mcsamplerGPU ")
    mcsampler_gpu_ok = False

import RIFT.likelihood.priors_utils as priors_utils
import RIFT.misc.xmlutils as xmlutils

#from lalinference.bayestar import fits as bfits

__author__ = "Evan Ochsner <evano@gravity.phys.uwm.edu>, Chris Pankow <pankow@gravity.phys.uwm.edu>, R. O'Shaughnessy"

#
# Pinnable parameters -- for command line processing
#
LIKELIHOOD_PINNABLE_PARAMS = ["right_ascension", "declination", "psi", "distance", "phi_orb", "t_ref", "inclination"]

def get_pinned_params(opts):
    """
    Retrieve a dictionary of user pinned parameters and their pin values.
    """
    return dict([(p,v) for p, v in opts.__dict__.items() if p in LIKELIHOOD_PINNABLE_PARAMS and v is not None]) 

def get_unpinned_params(opts, params):
    """
    Retrieve a set of unpinned parameters.
    """
    return params - set([p for p, v in opts.__dict__.items() if p in LIKELIHOOD_PINNABLE_PARAMS and v is not None])

#
# Option parsing
#

optp = OptionParser()
optp.add_option("-c", "--cache-file", default=None, help="LIGO cache file containing all data needed.")
optp.add_option("-C", "--channel-name", action="append", help="instrument=channel-name, e.g. H1=FAKE-STRAIN. Can be given multiple times for different instruments.")
optp.add_option("-p", "--psd-file", action="append", help="instrument=psd-file, e.g. H1=H1_PSD.xml.gz. Can be given multiple times for different instruments.")
optp.add_option("-k", "--skymap-file", help="Use skymap stored in given FITS file.")
optp.add_option("-x", "--coinc-xml", help="gstlal_inspiral XML file containing coincidence information.")
optp.add_option("-I", "--sim-xml", help="XML file containing mass grid to be evaluated")
optp.add_option("-E", "--event", default=0,type=int, help="Event number used for this run")
optp.add_option("--n-events-to-analyze", default=1,type=int, help="Number of events to analyze from this XML")
optp.add_option("--soft-fail-event-range",action='store_true',help='Soft failure (exit 0) if event ID is out of range. This happens in pipelines, if we have pre-built a DAG attempting to analyze more points than we really have')
optp.add_option("-f", "--reference-freq", type=float, default=100.0, help="Waveform reference frequency. Required, default is 100 Hz.")
optp.add_option("--fmin-template", dest='fmin_template', type=float, default=40, help="Waveform starting frequency.  Default is 40 Hz. Also equal to starting frequency for integration") 
optp.add_option("--fmin-template-correct-for-lmax",action='store_true',help="Modify amount of data selected, waveform starting frequency to account for l-max, to better insure all requested modes start within the targeted band")
optp.add_option("--fmin-ifo", action='append' , help="Minimum frequency for each IFO. Implemented by setting the PSD=0 below this cutoff. Use with care.") 
#optp.add_option("--nr-params",default=None, help="List of specific NR parameters and groups (and masses?) to use for the grid.")
#optp.add_option("--nr-index",type=int,default=-1,help="Index of specific NR simulation to use [integer]. Mass used: mtot= m1+m2")
optp.add_option('--nr-group', default=None,help="If using a *ssingle specific simulation* specified on the command line, provide it here")
optp.add_option('--nr-param', default=None,help="If using a *ssingle specific simulation* specified on the command line, provide it here")
optp.add_option("--nr-lookup",action='store_true', help=" Look up parameters from an NR catalog, instead of using the approximant specified")
optp.add_option("--nr-lookup-group",action='append', help="Restriction on 'group' for NR lookup")
optp.add_option("--nr-hybrid-use",action='store_true',help="Enable use of NR (or ROM!) hybrid, using --approx as the default approximant and with a frequency fmin")
optp.add_option("--nr-hybrid-method",default="taper_add",help="Hybridization method for NR (or ROM!).  Passed through to LALHybrid. pseudo_aligned_from22 will provide ad-hoc higher modes, if the early-time hybridization model only includes the 22 mode")
optp.add_option("--rom-group",default=None)
optp.add_option("--rom-param",default=None)
optp.add_option("--rom-use-basis",default=False,action='store_true',help="Use the ROM basis for inner products.")
optp.add_option("--rom-limit-basis-size-to",default=None,type=int)
optp.add_option("--rom-integrate-intrinsic",default=False,action='store_true',help='Integrate over intrinsic variables. REQUIRES rom_use_basis at present. ONLY integrates in mass ratio as present')
optp.add_option("--nr-perturbative-extraction",default=False,action='store_true')
optp.add_option("--nr-perturbative-extraction-full",default=False,action='store_true')
optp.add_option("--nr-use-provided-strain",default=False,action='store_true')
optp.add_option("--no-memory",default=False,action='store_true', help="At present, turns off m=0 modes. Use with EXTREME caution only if requested by model developer")
optp.add_option("--restricted-mode-list-file",default=None,help="A list of ALL modes to use in likelihood. Incredibly dangerous. Only use when comparing with models which provide restricted mode sets, or otherwise to isolate the effect of subsets of modes on the whole")
optp.add_option("--use-external-EOB",default=False,action='store_true')
optp.add_option("--maximize-only",default=False, action='store_true',help="After integrating, attempts to find the single best fitting point")
optp.add_option("--dump-lnL-time-series",default=False, action='store_true',help="(requires --sim-xml) Dump lnL(t) at the injected parameters")
optp.add_option("-a", "--approximant", default="TaylorT4", help="Waveform family to use for templates. Any approximant implemented in LALSimulation is valid.")
optp.add_option("-A", "--amp-order", type=int, default=0, help="Include amplitude corrections in template waveforms up to this e.g. (e.g. 5 <==> 2.5PN), default is Newtonian order.")
optp.add_option("--l-max", type=int, default=2, help="Include all (l,m) modes with l less than or equal to this value.")
optp.add_option("-s", "--data-start-time", type=float, default=None, help="GPS start time of data segment. If given, must also give --data-end-time. If not given, sane start and end time will automatically be chosen.")
optp.add_option("-e", "--data-end-time", type=float, default=None, help="GPS end time of data segment. If given, must also give --data-start-time. If not given, sane start and end time will automatically be chosen.")
optp.add_option("--data-integration-window-half",default=75*1e-3,type=float,help="Only change this window size if you are an expert. The window for time integration is -/+ this quantity around the event time")
optp.add_option("-F", "--fmax", type=float, help="Upper frequency of signal integration. Default is use PSD's maximum frequency.")
optp.add_option("--srate",default=16384,type=int,help="Sampling rate. Change ONLY IF YOU ARE ABSOLUTELY SURE YOU KNOW WHAT YOU ARE DOING.")
optp.add_option("-t", "--event-time", type=float, help="GPS time of the event --- probably the end time. Required if --coinc-xml not given.")
optp.add_option("-i", "--inv-spec-trunc-time", type=float, default=8., help="Timescale of inverse spectrum truncation in seconds (Default is 8 - give 0 for no truncation)")
optp.add_option("-w", "--window-shape", type=float, default=0, help="Shape of Tukey window to apply to data (default is no windowing)")
optp.add_option("--psd-window-shape", type=float, default=0, help="Shape of Tukey window that *was* applied to the PSD being passed. If nonzero, we will rescale the PSD by the ratio of window shape results")
optp.add_option("-m", "--time-marginalization", action="store_true", help="Perform marginalization over time via direct numerical integration. Default is false.")
optp.add_option("--vectorized", action="store_true", help="Perform manipulations of lm and timeseries using numpy arrays, not LAL data structures.  (Combine with --gpu to enable GPU use, where available)")
optp.add_option("--gpu", action="store_true", help="Perform manipulations of lm and timeseries using numpy arrays, CONVERTING TO GPU when available. You MUST use this option with --vectorized (otherwise it is a no-op). You MUST have a suitable version of cupy installed, your cuda operational, etc")
optp.add_option("--force-gpu-only", action="store_true", help="Hard fail if no GPU present (assessed by cupy not loading)")
optp.add_option("--force-xpy", action="store_true", help="Use the xpy code path.  Use with --vectorized --gpu to use the fallback CPU-based code path. Useful for debugging.")
optp.add_option("-o", "--output-file", help="Save result to this file.")
optp.add_option("-O", "--output-format", default='xml', help="[xml|hdf5]")
optp.add_option("-S", "--save-samples", action="store_true", help="Save sample points to output-file. Requires --output-file to be defined.")
optp.add_option("-L", "--save-deltalnL", type=float, default=float("Inf"), help="Threshold on deltalnL for points preserved in output file.  Requires --output-file to be defined")
optp.add_option("-P", "--save-P", type=float,default=0.1, help="Threshold on cumulative probability for points preserved in output file.  Requires --output-file to be defined")
optp.add_option("--internal-hard-fail-on-error",action='store_true',help='If true, fails with exit code 1 if any point is unsuccessful')
optp.add_option("--internal-make-empty-file-on-error",action='store_true',help='If true, failed points generate empty output file. Protects against OSG workflow problems')
optp.add_option("--verbose",action='store_true')
optp.add_option("--save-eccentricity", action="store_true")
#
# Add the integration options
#
integration_params = OptionGroup(optp, "Integration Parameters", "Control the integration with these options.")
# Default is actually None, but that tells the integrator to go forever or until n_eff is hit.
integration_params.add_option("--n-max", type=int, help="Total number of samples points to draw. If this number is hit before n_eff, then the integration will terminate. Default is 'infinite'.",default=1e7)
integration_params.add_option("--n-eff", type=int, default=100, help="Total number of effective samples points to calculate before the integration will terminate. Default is 100")
integration_params.add_option("--fairdraw-extrinsic-output", action='store_true' , help="Output is fair draw, rather than being comprehensive")
integration_params.add_option("--n-chunk", type=int, help="Chunk'.",default=10000)
integration_params.add_option("--convergence-tests-on",default=False,action='store_true')
integration_params.add_option("--seed", type=int, help="Random seed to use. Default is to not seed the RNG.")
integration_params.add_option("--no-adapt", action="store_true", help="Turn off adaptive sampling. Adaptive sampling is on by default.")
integration_params.add_option("--no-adapt-distance", action="store_true", help="Turn off adaptive sampling, just for distance. Adaptive sampling is on by default.")
integration_params.add_option("--no-adapt-after-first",action='store_true',help="Disables adaptation after first iteration with significant lnL")
integration_params.add_option("--adapt-weight-exponent", type=float, default=1.0, help="Exponent to use with weights (likelihood integrand) when doing adaptive sampling. Used in tandem with --adapt-floor-level to prevent overconvergence. Default is 1.0.")
integration_params.add_option("--adapt-floor-level", type=float, default=0.1, help="Floor to use with weights (likelihood integrand) when doing adaptive sampling. This is necessary to ensure the *sampling* prior is non zero during adaptive sampling and to prevent overconvergence. Default is 0.1 (no floor)")
integration_params.add_option("--adapt-adapt",action='store_true',help="Adapt the tempering exponent")
integration_params.add_option("--adapt-log",action='store_true',help="Use a logarithmic tempering exponent")
integration_params.add_option("--interpolate-time", default=False,help="If using time marginalization, compute using a continuously-interpolated array. (Default=false)")
integration_params.add_option("--d-prior",default='Euclidean' ,type=str,help="Distance prior for dL.  Options are dL^2 (Euclidean) and 'pseudo-cosmo'  .")
integration_params.add_option("--d-max", default=10000,type=float,help="Maximum distance in volume integral. Used to SET THE PRIOR; changing this value changes the numerical answer.")
integration_params.add_option("--d-min", default=1,type=float,help="Minimum distance in volume integral. Used to SET THE PRIOR; changing this value changes the numerical answer.")
integration_params.add_option("--declination-cosine-sampler",action='store_true',help="If specified, the parameter used for declination is cos(dec), not dec")
integration_params.add_option("--inclination-cosine-sampler",action='store_true',help="If specified, the parameter used for inclination is cos(dec), not dec")
integration_params.add_option("--manual-logarithm-offset",type=float,default=0,help="Target value of logarithm lnL. Integrand is reduced by exp(-manual_logarithm_offset).  Important for high-SNR sources!   Should be set dynamically")
integration_params.add_option("--sampler-method",default="adaptive_cartesian",help="adaptive_cartesian|GMM|adaptive_cartesian_gpu")
integration_params.add_option("--sampler-xpy",default=None,help="numpy|cupy  if the adaptive_cartesian_gpu sampler is active, use that.")
optp.add_option_group(integration_params)

#
# Add the intrinsic parameters
#
intrinsic_params = OptionGroup(optp, "Intrinsic Parameters", "Intrinsic parameters (e.g component mass) to use.")
#intrinsic_params.add_option("--pin-distance-to-sim",action='store_true', help="Pin *distance* value to sim entry. Used to enable source frame reconstruction with NR.")
intrinsic_params.add_option("--mass1", type=float, help="Value of first component mass, in solar masses. Required if not providing coinc tables.")
intrinsic_params.add_option("--mass2", type=float, help="Value of second component mass, in solar masses. Required if not providing coinc tables.")
intrinsic_params.add_option("--eff-lambda", type=float, help="Value of effective tidal parameter. Optional, ignored if not given.")
intrinsic_params.add_option("--deff-lambda", type=float, help="Value of second effective tidal parameter. Optional, ignored if not given")
optp.add_option_group(intrinsic_params)


#
# Add options to integrate over intrinsic parameters.  Same conventions as util_ManualOverlapGrid.py.  
# Parameters have special names, and we adopt priors that use those names.
# NOTE: Only 'q' implemented
#
intrinsic_int_params = OptionGroup(optp, "Intrinsic integrated parameters", "Intrinsic parameters to integrate over. ONLY currently used with ROM version")
intrinsic_int_params.add_option("--parameter",action='append')
intrinsic_int_params.add_option("--parameter-range",action='append',type=str)
intrinsic_int_params.add_option("--adapt-intrinsic",action='store_true')
optp.add_option_group(intrinsic_int_params)


#
# Add the pinnable parameters
#
pinnable = OptionGroup(optp, "Pinnable Parameters", "Specifying these command line options will pin the value of that parameter to the specified value with a probability of unity.")
for pin_param in LIKELIHOOD_PINNABLE_PARAMS:
    option = "--" + pin_param.replace("_", "-")
    pinnable.add_option(option, type=float, help="Pin the value of %s." % pin_param)
optp.add_option_group(pinnable)

opts, args = optp.parse_args()

if opts.gpu and opts.force_gpu_only and not cupy_success:
  print("GPU requested but no GPU/cupy found: Hard fail by request")
  sys.exit(35)  # unique code to make resuming due to this error easier to identify

if opts.gpu and xpy_default is numpy:
    print( " Override --gpu  (not available);  use --force-xpy to require the identical code path is used (with xpy =[np|cupy]")
    opts.gpu=False
    if opts.force_xpy:
        opts.gpu=True

manual_avoid_overflow_logarithm=opts.manual_logarithm_offset

deltaT = None
fSample= opts.srate # change sampling rate
if not(fSample is None):
  deltaT =1./fSample


# Load in restricted mode set, if available
restricted_mode_list=None
if not(opts.restricted_mode_list_file is None):
    modes =numpy.loadtxt(opts.restricted_mode_list_file,dtype=int) # columns are l m.  Must contain all. Only integers obviously
    restricted_mode_list = [ (l,m) for l,m in modes]
    print( " RESTRICTED MODE LIST target :", restricted_mode_list)

intrinsic_param_names = opts.parameter
valid_intrinsic_param_names = ['q']
if intrinsic_param_names:
 for param in intrinsic_param_names:
    
    # Check if in the valid list
    if not(param in valid_intrinsic_param_names):
            print( ' Invalid param ', param, ' not in ', valid_intrinsic_param_names)
            sys.exit(1)
    param_ranges = []
    if len(intrinsic_param_names) == len(opts.parameter_range):
        param_ranges = numpy.array(map(eval, opts.parameter_range))
        # Rescale mass-dependent ranges to SI units
        for p in ['mc', 'm1', 'm2', 'mtot']:
          if p in intrinsic_param_names:
            indx = intrinsic_param_names.index(p)
            param_ranges[indx]= numpy.array(param_ranges[indx])* lal.MSUN_SI



# Check both or neither of --data-start/end-time given
if opts.data_start_time is None and opts.data_end_time is not None:
    raise ValueError("You must provide both or neither of --data-start-time and --data-end-time.")
if opts.data_end_time is None and opts.data_start_time is not None:
    raise ValueError("You must provide both or neither of --data-start-time and --data-end-time.")

#
# Import NR grid
#
NR_template_group=None
NR_template_param=None
if opts.nr_group and opts.nr_param:
    import NRWaveformCatalogManager3 as nrwf
    NR_template_group = opts.nr_group
    if nrwf.internal_ParametersAreExpressions[NR_template_group]:
        NR_template_param = eval(opts.nr_param)
    else:
        NR_template_param = opts.nr_param


#
# Hardcoded variables
#
template_min_freq = opts.fmin_template # minimum frequency of template
#t_ref_wind = 50e-3 # Interpolate in a window +/- this width about event time. 
t_ref_wind = opts.data_integration_window_half
T_safety = 2. # Safety buffer (in sec) for wraparound corruption

#
# Inverse spectrum truncation control
#
T_spec = opts.inv_spec_trunc_time
if T_spec == 0.: # Do not do inverse spectrum truncation
    inv_spec_trunc_Q = False
    T_safety += 8. # Add a bit more safety buffer in this case
else:
    inv_spec_trunc_Q = True

#
# Integrator options
#
n_max = opts.n_max # Max number of extrinsic points to evaluate at
n_eff = opts.n_eff # Effective number of points evaluated

#
# Initialize the RNG, if needed
#
# TODO: Do we seed a given instance of the integrator, or set it for all
# or both?
if opts.seed is not None:
    numpy.random.seed(opts.seed)


if opts.event_time is not None:
    event_time = glue.lal.LIGOTimeGPS(opts.event_time)
    print( "Event time from command line: %s" % str(event_time))
else:
    print( " Error! ")
    sys.exit(1)

#
# Template descriptors
#

fiducial_epoch = lal.LIGOTimeGPS()
fiducial_epoch = event_time.seconds + 1e-9*event_time.nanoseconds   # no more direct access to gpsSeconds

# Struct to hold template parameters
P_list = None
if opts.sim_xml:
    print( "====Loading injection XML:", opts.sim_xml, opts.event, " =======")
    P_list = lalsimutils.xml_to_ChooseWaveformParams_array(str(opts.sim_xml))
    if  len(P_list) < opts.event+opts.n_events_to_analyze:
        print( " Event list of range; soft exit")
        sys.exit(0)
    P_list = P_list[opts.event:(opts.event+opts.n_events_to_analyze)]
    for P in P_list:
        P.radec =False  # do NOT propagate the epoch later
        P.fref = opts.reference_freq
        P.fmin = template_min_freq
        P.tref = fiducial_epoch  # the XML table
        m1 = P.m1/lal.MSUN_SI
        m2 =P.m2/lal.MSUN_SI
        lambda1, lambda2 = P.lambda1,P.lambda2
        P.dist = factored_likelihood.distMpcRef * 1.e6 * lal.PC_SI   # use *nonstandard* distance
        P.phi=0.0
        P.psi=0.0
        P.incl = 0.0       # only works for aligned spins. Be careful.
        P.fref = opts.reference_freq
        if opts.approximant != "TaylorT4": # not default setting
            P.approx = lalsimutils.lalsim.GetApproximantFromString(opts.approximant)  # allow user to override the approx setting. Important for NR followup, where no approx set in sim_xml!
        if opts.approximant == "EccentricTD":
            P.phaseO = 3

    P = P_list[0]  # Load in the physical parameters of the injection.  


# User requested bounds for data segment
if not (opts.data_start_time == None) and  not (opts.data_end_time == None):
    start_time =  opts.data_start_time
    end_time =  opts.data_end_time
    print( "Fetching data segment with start=", start_time)
    print( "                             end=", end_time)

# Automatically choose data segment bounds so region of interest isn't corrupted
# FIXME: Use estimate, instead of painful full waveform generation call here.
else:
    approxTmp = P.approx
    LmaxEff = 2
    if opts.fmin_template_correct_for_lmax:
      LmaxEff=opts.l_max
    T_tmplt = lalsimutils.estimateWaveformDuration(P,LmaxEff) + 4  # Much more robust than the previous (slow, prone-to-crash approach)
#     if (m1+m2)<30:
#         P.approx = lalsimutils.lalsim.TaylorT4  # should not impact length much.  IMPORTANT because EOB calls will fail at the default sampling rate
#         htmplt = lalsimutils.hoft(P)   # Horribly wasteful waveform gneeration solely to estimate duration.  Will also crash if spins are used.
#         P.approx = approxTmp
#         T_tmplt = - float(htmplt.epoch)    # ASSUMES time returned by hlm is a RELATIVE time that does not add tref. P.radec=False !
#         print fiducial_epoch, event_time, htmplt.epoch
#     else:
#         print " Using estimate for waveform length in a high mass regime: beware!"
#         T_tmplt = lalsimutils.estimateWaveformDuration(P) + 4
    T_seg = T_tmplt + T_spec + T_safety # Amount before and after event time
#    if opts.use_external_EOB:
#        T_seg *=2   # extra safety factor
    start_time = float(event_time) - T_seg
    end_time = float(event_time) + T_seg
    print( "Fetching data segment with start=", start_time)
    print( "                             end=", end_time)
    print( "\t\tEvent time is: ", float(event_time))
    print( "\t\tT_seg is: ", T_seg)

#
# Load in data and PSDs
#
data_dict, psd_dict = {}, {}

for inst, chan in map(lambda c: c.split("="), opts.channel_name):
    print( "Reading channel %s from cache %s" % (inst+":"+chan, opts.cache_file))
    data_dict[inst] = lalsimutils.frame_data_to_non_herm_hoff(opts.cache_file,
            inst+":"+chan, start=start_time, stop=end_time,
            window_shape=opts.window_shape,deltaT=deltaT)
    print( "Frequency binning: %f, length %d" % (data_dict[inst].deltaF,
            data_dict[inst].data.length) )

flow_ifo_dict = {}
if opts.fmin_ifo:
 for inst, freq_str in map(lambda c: c.split("="), opts.fmin_ifo):
    freq_low_here = float(freq_str)
    print( "Reading low frequency cutoff for instrument %s from %s" % (inst, freq_str), freq_low_here)
    flow_ifo_dict[inst] = freq_low_here

for inst, psdf in map(lambda c: c.split("="), opts.psd_file):
    print( "Reading PSD for instrument %s from %s" % (inst, psdf))
    psd_dict[inst] = lalsimutils.get_psd_series_from_xmldoc(psdf, inst)

    deltaF = data_dict[inst].deltaF
    psd_dict[inst] = lalsimutils.resample_psd_series(psd_dict[inst], deltaF)
    print( "PSD deltaF after interpolation %f" % psd_dict[inst].deltaF)

    # Implement PSD window rescaling: see T1900249
    if opts.psd_window_shape > 0 or opts.window_shape > 0:
      # Note this is just a constant scale factor in the likelihood, but it *can* have some effect on parameters if one windowing size is small (as it stretches/compresses the likleihood)
      window_fac_psd = lalsimutils.psd_windowing_factor(opts.psd_window_shape, len(psd_dict[inst].data.data))  # assume the windowing factor IS accounted for, and we have to undo it
      window_fac_data = lalsimutils.psd_windowing_factor(opts.window_shape, len(data_dict[inst].data.data))
      psd_dict[inst].data.data *= window_fac_data/window_fac_psd   # scale the windowing factor to be appropriate for our actual data window: the PSD was measured using a windowing different than the one we use

    # implement cutoff.  
    if inst in flow_ifo_dict.keys():
        if isinstance(psd_dict[inst], lal.REAL8FrequencySeries):
            psd_fvals = psd_dict[inst].f0 + deltaF*numpy.arange(psd_dict[inst].data.length)
            psd_dict[inst].data.data[ psd_fvals < flow_ifo_dict[inst]] = 0 # 
        else:
            print( 'FAIL on PSD import')
            sys.exit(1)
#        elif isinstance(psd_dict[inst], pylal.xlal.datatypes.real8frequencyseries.REAL8FrequencySeries):  # for backward compatibility
#            psd_fvals = psd_dict[inst].f0 + deltaF*numpy.arange(len(psd_dict[inst].data))
#            psd_dict[inst].data[psd_fvals < ifo_dict[inst]] =0


    assert psd_dict[inst].deltaF == deltaF

    # Highest freq. at which PSD is defined
    # if isinstance(psd_dict[inst],
    #         pylal.xlal.datatypes.real8frequencyseries.REAL8FrequencySeries):
    #     fmax = psd_dict[inst].f0 + deltaF * (len(psd_dict[inst].data) - 1)
    if isinstance(psd_dict[inst], lal.REAL8FrequencySeries):
        fmax = psd_dict[inst].f0 + deltaF * (psd_dict[inst].data.length - 1)

    # Assert upper limit of IP integral does not go past where PSD defined
    assert opts.fmax is None or opts.fmax<= fmax
    # Allow us to target a smaller upper limit than provided by the PSD. Important for numerical PSDs that turn over at high frequency
    if opts.fmax and opts.fmax < fmax:
        fmax = opts.fmax # fmax is now the upper freq. of IP integral

# Ensure data and PSDs keyed to same detectors
if sorted(psd_dict.keys()) != sorted(data_dict.keys()):
    print >>sys.stderr, "Got a different set of instruments based on data and PSDs provided."

# Ensure waveform has same sample rate, padded length as data
#
# N.B. This assumes all detector data has same sample rate, length
#
# data_dict holds 2-sided FrequencySeries, so their length is the same as
# that of the original TimeSeries that was FFT'd = Nsamples
# Also, deltaF = 1/T, with T = the duration (in sec) of the original TimeSeries
# Therefore 1/(data.length*deltaF) = T/Nsamples = deltaT
P.deltaT = 1./ (data_dict[list(data_dict.keys())[0]].data.length * deltaF)
P.deltaF = deltaF
for Psig in P_list:
  Psig.deltaT = P.deltaT
  Psig.deltaf = P.deltaF



#
# Set up parameters and bounds
#

# PROBLEM: if too large, you can MISS a source. Does NOT need to be fixed for all masses *IF* the problem really has strong support
dmin = opts.d_min    # min distance
dmax = opts.d_max  # max distance FOR ANY SOURCE EVER. EUCLIDEAN



dmax_sampling_guess = dmax
distBoundGuess = dmax

print( "Recommended distance for sampling ", dmax_sampling_guess, " and probably near ", distBoundGuess, " smaller than  ", dmax)
print( "    (recommendation not yet used) ")
param_limits = { "psi": (0, 2*numpy.pi),
    "phi_orb": (0, 2*numpy.pi),
    "distance": (dmin, dmax),   # CAN LEAD TO CATASTROPHIC FAILURE if dmax is too large (adaptive fails - too few bins)
    "right_ascension": (0, 2*numpy.pi),
    "declination": (-numpy.pi/2, numpy.pi/2),
    "t_ref": (-t_ref_wind, t_ref_wind),
    "inclination": (0, numpy.pi)
}

#
# Parameter integral sampling strategy
#
params = {}
sampler = mcsampler.MCSampler()
if opts.sampler_method == "adaptive_cartesian_gpu":
    sampler = mcsamplerGPU.MCSampler()
    sampler.xpy = xpy_default
    sampler.identity_convert=identity_convert
    mcsampler  = mcsamplerGPU  # force use of routines in that file, for properly configured GPU-accelerated code as needed

    if opts.sampler_xpy == "numpy":
      mcsampler.set_xpy_to_numpy()
      sampler.xpy= numpy
      sampler.identity_convert= lambda x: x
if opts.sampler_method == "GMM":
    sampler = mcsamplerEnsemble.MCSampler()


#
# Psi -- polarization angle
# sampler: uniform in [0, pi)
#
psi_sampler = mcsampler.ret_uniform_samp_vector_alt( 
    param_limits["psi"][0], param_limits["psi"][1])
psi_sampler_cdf_inv = functools.partial(mcsampler.uniform_samp_cdf_inv_vector, 
    param_limits["psi"][0], param_limits["psi"][1])
sampler.add_parameter("psi", 
    pdf = psi_sampler, 
    cdf_inv = psi_sampler_cdf_inv, 
    left_limit = param_limits["psi"][0], 
    right_limit = param_limits["psi"][1],
    prior_pdf = mcsampler.uniform_samp_psi)

#
# Phi - orbital phase
# sampler: uniform in [0, 2*pi)
#
phi_sampler = mcsampler.ret_uniform_samp_vector_alt( 
    param_limits["phi_orb"][0], param_limits["phi_orb"][1])
phi_sampler_cdf_inv = functools.partial(mcsampler.uniform_samp_cdf_inv_vector, 
    param_limits["phi_orb"][0], param_limits["phi_orb"][1])
sampler.add_parameter("phi_orb",
    pdf = phi_sampler,
    cdf_inv = phi_sampler_cdf_inv,
    left_limit = param_limits["phi_orb"][0], 
    right_limit = param_limits["phi_orb"][1],
    prior_pdf = mcsampler.uniform_samp_phase)

#
# inclination - angle of system angular momentum with line of sight
# sampler: cos(incl) uniform in [-1, 1)
#

adapt_extra_extrinsic=False
if opts.sampler_method == "adaptive_cartesian_gpu":  # this is a better/more stable/faster adaptive code, trust it to adapt in more extrinsic dimensions
  adapt_extra_extrinsic=True

if not opts.inclination_cosine_sampler:
 incl_sampler = mcsampler.cos_samp_vector # this is NOT dec_samp_vector, because the angular zero point is different!
 incl_sampler_cdf_inv = mcsampler.cos_samp_cdf_inv_vector
 sampler.add_parameter("inclination", 
    pdf = incl_sampler, 
    cdf_inv = incl_sampler_cdf_inv, 
    left_limit = param_limits["inclination"][0], 
    right_limit = param_limits["inclination"][1],
    prior_pdf = mcsampler.uniform_samp_theta)  # do NOT adapt if we are using a prior that goes to zero at the edges
else:
 incl_sampler = mcsampler.ret_uniform_samp_vector_alt(-1.0,1.0)
 incl_sampler_cdf_inv = lambda x: x*2.0-1.  #functools.partial(mcsampler.uniform_samp_cdf_inv_vector,-1,1) 
 sampler.add_parameter("inclination", 
    pdf = incl_sampler, 
    cdf_inv = incl_sampler_cdf_inv, 
    left_limit = -1, 
    right_limit = 1,
    prior_pdf = incl_sampler,
    adaptive_sampling=adapt_extra_extrinsic)

#
# Distance - luminosity distance to source in parsecs
# sampler: uniform distance over [dmin, dmax), adaptive sampling
#
dist_sampler = mcsampler.ret_uniform_samp_vector_alt( param_limits["distance"][0], param_limits["distance"][1])
dist_sampler_cdf_inv = functools.partial(mcsampler.uniform_samp_cdf_inv_vector,     param_limits["distance"][0], param_limits["distance"][1])
#dist_sampler=functools.partial( mcsampler.uniform_samp_withfloor_vector, numpy.min([distBoundGuess,param_limits["distance"][1]]), param_limits["distance"][1], 0.001)
dist_prior_pdf =   lambda x: x**2/(param_limits["distance"][1]**3/3.                   - param_limits["distance"][0]**3/3.) 
if opts.d_prior == 'pseudo_cosmo':
  nm = priors_utils.dist_prior_pseudo_cosmo_eval_norm(param_limits["distance"][0],param_limits["distance"][1])
  dist_prior_pdf =functools.partial( priors_utils.dist_prior_pseudo_cosmo, nm=nm)
#dist_sampler_cdf_inv=None
sampler.add_parameter("distance", 
    pdf = dist_sampler, 
    cdf_inv = dist_sampler_cdf_inv,
    left_limit = param_limits["distance"][0], 
    right_limit = param_limits["distance"][1],
    prior_pdf = dist_prior_pdf,
    adaptive_sampling = adapt_extra_extrinsic and not (opts.no_adapt or opts.no_adapt_distance))


#
# Intrinsic parameters
#
sampler_lookup = {}
sampler_inv_lookup = {}
sampler_lookup['q'] = mcsampler.q_samp_vector   # only one intrinsic parameter possible
sampler_lookup['M'] = mcsampler.M_samp_vector
sampler_inv_lookup['q'] = mcsampler.q_cdf_inv_vector
sampler_inv_lookup['M'] = None # mcsampler.M_cdf_inv_vector

if opts.rom_use_basis and opts.rom_integrate_intrinsic:
 for p in intrinsic_param_names:
    indx = intrinsic_param_names.index(p)
    qmin,qmax = param_ranges[indx]
    q_pdf =mcsampler.ret_uniform_samp_vector_alt( qmin, qmax)   # sample uniformly by default
#    q_pdf =functools.partial(sampler_lookup[p], qmin, qmax)   # sample uniformly by default
    q_cdf_inv = functools.partial(mcsampler.uniform_samp_cdf_inv_vector, qmin,qmax) # sample uniformly by default
#    q_cdf_inv= None
    q_pdf_prior = functools.partial(sampler_lookup[param], qmin, qmax)  # true prior, from lookup
    sampler.add_parameter(p, 
        pdf=q_pdf, 
        cdf_inv=q_cdf_inv, 
        left_limit=qmin, right_limit=qmax,prior_pdf=q_pdf_prior, 
        adaptive_sampling = opts.adapt_intrinsic)



if opts.skymap_file is not None:
    #
    # Right ascension and declination -- use a provided skymap
    #
    smap, _ = bfits.read_sky_map(opts.skymap_file)
    # FIXME: Uncomment for 'mixed' map
    #smap = 0.9*smap + 0.1*numpy.ones(len(smap))/len(smap)
    ss_sampler = mcsampler.HealPixSampler(smap)
    #isotropic_2d_sampler = numpy.vectorize(lambda dec, ra: mcsampler.dec_samp_vector(dec)/2/numpy.pi)
    isotropic_bstar_sampler = numpy.vectorize(lambda dec, ra: 1.0/len(smap))

    # FIXME: Should the left and right limits be modified?
    sampler.add_parameter(("declination", "right_ascension"), 
        pdf = ss_sampler.pseudo_pdf,
        cdf_inv = ss_sampler.pseudo_cdf_inverse, 
        left_limit = (param_limits["declination"][0], param_limits["right_ascension"][0]),
        right_limit = (param_limits["declination"][1], param_limits["right_ascension"][1]),
        prior_pdf = isotropic_bstar_sampler)

else:
    #
    # Right ascension - angle in radians from prime meridian plus hour angle
    # sampler: uniform in [0, 2pi), adaptive sampling
    #
    ra_sampler = mcsampler.ret_uniform_samp_vector_alt(
        param_limits["right_ascension"][0], param_limits["right_ascension"][1])
    ra_sampler_cdf_inv = functools.partial(mcsampler.uniform_samp_cdf_inv_vector,
        param_limits["right_ascension"][0], param_limits["right_ascension"][1])
    sampler.add_parameter("right_ascension", 
        pdf = ra_sampler, 
        cdf_inv = ra_sampler_cdf_inv, 
        left_limit = param_limits["right_ascension"][0],
        right_limit =  param_limits["right_ascension"][1],
        prior_pdf = mcsampler.uniform_samp_phase,
        adaptive_sampling = not opts.no_adapt)

    #
    # declination - angle in radians from the north pole piercing the celestial
    # sky sampler: cos(dec) uniform in [-1, 1), adaptive sampling
    #
    if not opts.declination_cosine_sampler:
     dec_sampler = mcsampler.dec_samp_vector
     dec_sampler_cdf_inv = mcsampler.dec_samp_cdf_inv_vector
     sampler.add_parameter("declination", 
        pdf = dec_sampler, 
        cdf_inv = dec_sampler_cdf_inv, 
        left_limit = param_limits["declination"][0], 
        right_limit = param_limits["declination"][1],
        prior_pdf = mcsampler.uniform_samp_dec,
        adaptive_sampling = not opts.no_adapt)
    else:
     # Sample uniformly in cos(polar_theta), =1 for north pole, -1 for south pole.
     # Propagate carefully in conversions: time of flight libraries use RA,DEC
     dec_sampler = mcsampler.ret_uniform_samp_vector_alt(-1.0,1.0)
     dec_sampler_cdf_inv = lambda x: x*2.0-1. # functools.partial(mcsampler.uniform_samp_cdf_inv_vector,-1,1)
     sampler.add_parameter("declination", 
        pdf = dec_sampler, 
        cdf_inv = dec_sampler_cdf_inv, 
        left_limit = -1, 
        right_limit = 1,
        prior_pdf = dec_sampler,
        adaptive_sampling = not opts.no_adapt)

if not opts.time_marginalization:
    #
    # tref - GPS time of geocentric end time
    # sampler: uniform in +/-2 ms window around estimated end time 
    #
    tref_sampler = mcsampler.ret_uniform_samp_vector_alt(
        param_limits["t_ref"][0], param_limits["t_ref"][1])

    tref_sampler_cdf_inv = functools.partial(mcsampler.uniform_samp_cdf_inv_vector, 
                                            param_limits["t_ref"][0], param_limits["t_ref"][1])
    sampler.add_parameter("t_ref", 
                          pdf = tref_sampler, 
                          cdf_inv = None, 
                          left_limit = param_limits["t_ref"][0], 
                          right_limit = param_limits["t_ref"][1],
                          prior_pdf = functools.partial(mcsampler.uniform_samp_vector, param_limits["t_ref"][0], param_limits["t_ref"][1]))


#
# Determine pinned and non-pinned parameters
#

pinned_params = get_pinned_params(opts)
unpinned_params = get_unpinned_params(opts, sampler.params)
print( "{0:<25s} {1:>5s} {2:>5s} {3:>20s} {4:<10s}".format("parameter", "lower limit", "upper limit", "pinned?", "pin value"))
plen = len(sorted(sampler.params, key=lambda p: len(p))[-1])
for p in sampler.params:
    if p in pinned_params:
        pinned, value = True, "%1.3g" % pinned_params[p]
    else:
        pinned, value = False, ""

    if isinstance(p, tuple):
        for subp, subl, subr in zip(p, sampler.llim[p], sampler.rlim[p]):
            subp = subp + " "*min(0, plen-len(subp))
            print( "|{0:<25s} {1:>1.3g}   {2:>1.3g} {3:>20s} {4:<10s}".format(subp, subl, subr, str(False), ""))
    else:
        p = p + " "*min(0, plen-len(p))
        print( "{0:<25s} {1:>1.3g}   {2:>1.3g} {3:>20s} {4:<10s}".format(p, sampler.llim[p], sampler.rlim[p], str(pinned), value))

# Special case: t_ref is assumed to be relative to the epoch
if "t_ref" in pinned_params:
    pinned_params["t_ref"] -= float(fiducial_epoch)

#
# Provide convergence tests
# FIXME: Currently using hardcoded thresholds, poorly hand-tuned
#
test_converged = {}

#
# Merge options into one big ol' kwargs dict
#

pinned_params.update({ 
    # Iteration settings and termination conditions
    "n": min(opts.n_chunk, n_max), # Number of samples in a chunk
    "nmax": n_max, # Total number of samples to draw before termination
    "neff": n_eff, # Total number of effective samples to collect before termination

    "convergence_tests" : test_converged,    # Dictionary of convergence tests

    # Adaptive sampling settings
    "tempering_exp": opts.adapt_weight_exponent if not opts.no_adapt else 0.0, # Weights will be raised to this power to prevent overconvergence
    "tempering_log": opts.adapt_log,
    "tempering_adapt": opts.adapt_adapt,

    "floor_level": opts.adapt_floor_level if not opts.no_adapt else 0.0, # The new sampling distribution at the end of each chunk will be floor_level-weighted average of a uniform distribution and the (L^tempering_exp p/p_s)-weighted histogram of sampled points.
    "history_mult": 10, # Multiplier on 'n' - number of samples to estimate marginalized 1-D histograms
    "n_adapt": 100 if not opts.no_adapt else 0, # Number of chunks to allow adaption over

    # Verbosity settings
    "verbose": True, #not opts.rom_integrate_intrinsic, 
    "extremely_verbose": False, 

    # Sample caching
    "save_intg": opts.save_samples, # Cache the samples (and integrand values)?
    "igrand_threshold_deltalnL": opts.save_deltalnL, # Threshold on distance from max L to save sample
    "igrand_threshold_p": opts.save_P, # Threshold on cumulative probability contribution to cache sample
    "igrand_fairdraw_samples": opts.fairdraw_extrinsic_output,
    "igrand_fairdraw_samples_max": opts.n_eff
})
if opts.sampler_method == "adaptive_cartesian_gpu":
  pinned_params.update({"save_no_samples":True})   # do not exhaust GPU memory with MC samples!  
if opts.sampler_method == "GMM":
    n_step =pinned_params["n"]
    n_max_blocks = ((1.0*int(opts.n_max))/n_step)
    # pairing coordinates for adaptive integration: see definition of order below
    #    (distance,inclination)
    #    (ra, dec)
    #    (psi) (orb_phase)   # only do 1d adaptivity there, 
#    gmm_dict = {tuple([2]):None,tuple([4]):None,(3,5):None,(0,1):None} 
#    comp_dict = {tuple([2]):1,tuple([4]):1,(3,5):3,(0,1):4} 
    gmm_dict = {tuple([2]):None,tuple([4]):None,(3,5):None,tuple([0]):None,tuple([1]):None} 
    n_sky = 4
#    if len(psd_dict.keys()) > 2:
#        n_sky=2  # presumably we have localized better, avoid failure modes
    comp_dict = {tuple([2]):1,tuple([4]):1,(3,5):3,tuple([0]):n_sky,tuple([1]):n_sky} 

    extra_args = {'n_comp':comp_dict,'max_iter':n_max_blocks,'gmm_dict':gmm_dict}  # made up for now, should adjust
    pinned_params.update(extra_args)

    # cannot pin params at present
    # unpinned params MUST be full call signature of likelihood, IN ORDER
    if not(opts.time_marginalization):
      unpinned_params =['right_ascension','declination', 't_ref','phi_orb','inclination','psi','distance']
    else:
      unpinned_params =['right_ascension','declination', 'phi_orb','inclination','psi','distance']

def analyze_event(P_list, indx_event, data_dict, psd_dict, fmax, opts,inv_spec_trunc_Q=inv_spec_trunc_Q, T_spec=T_spec):
    nEvals=0
    P = P_list[indx_event]
    # Precompute
    t_window = 0.15
    rholms_intp, cross_terms, cross_terms_V,  rholms, rest=factored_likelihood.PrecomputeLikelihoodTerms(
            fiducial_epoch, t_window, P, data_dict, psd_dict, opts.l_max, fmax,
            False, inv_spec_trunc_Q, T_spec,
            NR_group=NR_template_group,NR_param=NR_template_param,
            use_external_EOB=opts.use_external_EOB,nr_lookup=opts.nr_lookup,nr_lookup_valid_groups=opts.nr_lookup_group,perturbative_extraction=opts.nr_perturbative_extraction,perturbative_extraction_full=opts.nr_perturbative_extraction_full,use_provided_strain=opts.nr_use_provided_strain,hybrid_use=opts.nr_hybrid_use,hybrid_method=opts.nr_hybrid_method,ROM_group=opts.rom_group,ROM_param=opts.rom_param,ROM_use_basis=opts.rom_use_basis,verbose=opts.verbose,quiet=not opts.verbose,ROM_limit_basis_size=opts.rom_limit_basis_size_to,no_memory=opts.no_memory,skip_interpolation=opts.vectorized)

    if opts.vectorized:
        lookupNKDict = {}
        lookupKNDict={}
        lookupKNconjDict={}
        ctUArrayDict = {}
        ctVArrayDict={}
        rholmArrayDict={}
        rholms_intpArrayDict={}
        epochDict={}
        for det in rholms_intp.keys():
            print( " Packing ", det)
            lookupNKDict[det],lookupKNDict[det], lookupKNconjDict[det], ctUArrayDict[det], ctVArrayDict[det], rholmArrayDict[det], rholms_intpArrayDict[det], epochDict[det] = factored_likelihood.PackLikelihoodDataStructuresAsArrays( rholms[det].keys(), rholms_intp[det], rholms[det], cross_terms[det],cross_terms_V[det])
            if opts.gpu and (not xpy_default is np):
                lookupNKDict[det] = cupy.asarray(lookupNKDict[det])
                rholmArrayDict[det] = cupy.asarray(rholmArrayDict[det])
                ctUArrayDict[det] = cupy.asarray(ctUArrayDict[det])
                ctVArrayDict[det] = cupy.asarray(ctVArrayDict[det])
                epochDict[det] = cupy.asarray(epochDict[det])



    # Likelihood
    if not opts.time_marginalization:

        def likelihood_function(right_ascension, declination, t_ref, phi_orb,
                inclination, psi, distance):

            dec = numpy.copy(declination).astype(numpy.float64)
            if opts.declination_cosine_sampler:
                dec = numpy.pi/2 - numpy.arccos(dec)
            incl = numpy.copy(inclination).astype(numpy.float64)
            if opts.inclination_cosine_sampler:
                incl = numpy.arccos(incl)

            # use EXTREMELY many bits
            lnL = numpy.zeros(right_ascension.shape,dtype=numpy.float128)
            i = 0
            for ph, th, tr, phr, ic, ps, di in zip(right_ascension, dec,
                    t_ref, phi_orb, incl, psi, distance):
                P.phi = ph # right ascension
                P.theta = th # declination
                P.tref = fiducial_epoch + tr # ref. time (rel to epoch for data taking)
                P.phiref = phr # ref. orbital phase
                P.incl = ic # inclination
                P.psi = ps # polarization angle
                P.dist = di* 1.e6 * lalsimutils.lsu_PC # luminosity distance

                lnL[i] = factored_likelihood.FactoredLogLikelihood(
                        P, rholms_intp, cross_terms, cross_terms_V, opts.l_max)
                i+=1

            return numpy.exp(lnL - manual_avoid_overflow_logarithm)

    else: # Sum over time for every point in other extrinsic params
     if not (opts.rom_integrate_intrinsic or opts.vectorized):
        def likelihood_function(right_ascension, declination, phi_orb, inclination,
                psi, distance):
            dec = numpy.copy(declination).astype(numpy.float64)  # get rid of 'object', and allocate space
            if opts.declination_cosine_sampler:
                dec = numpy.pi/2 - numpy.arccos(dec)
            incl = numpy.copy(inclination).astype(numpy.float64)
            if opts.inclination_cosine_sampler:
                incl = numpy.arccos(incl)

            # use EXTREMELY many bits
            lnL = numpy.zeros(right_ascension.shape,dtype=numpy.float128)
            i = 0
            tvals = numpy.linspace(-t_ref_wind,t_ref_wind,int((t_ref_wind)*2/P.deltaT))  # choose an array at the target sampling rate. P is inherited globally

            for ph, th, phr, ic, ps, di in zip(right_ascension, dec,
                    phi_orb, inclination, psi, distance):
                P.phi = ph # right ascension
                P.theta = th # declination
                P.tref = fiducial_epoch  # see 'tvals', above
                P.phiref = phr # ref. orbital phase
                P.incl = ic # inclination
                P.psi = ps # polarization angle
                P.dist = di* 1.e6 * lalsimutils.lsu_PC # luminosity distance


                lnL[i] = factored_likelihood.FactoredLogLikelihoodTimeMarginalized(tvals,
                        P, rholms_intp, rholms, cross_terms, cross_terms_V,                   
                        opts.l_max,interpolate=opts.interpolate_time)
                i+=1

            return numpy.exp(lnL -manual_avoid_overflow_logarithm)

     elif opts.vectorized: # use array-based multiplications, fewer for loops
        if (not opts.gpu):
          def likelihood_function(right_ascension, declination, phi_orb, inclination,
                psi, distance):
#            global nEvals
            tvals = numpy.linspace(-t_ref_wind,t_ref_wind,int((t_ref_wind)*2/P.deltaT))  # choose an array at the target sampling rate. P is inherited globally
            dec = numpy.copy(declination).astype(numpy.float64)
            if opts.declination_cosine_sampler:
              dec = numpy.pi/2 - numpy.arccos(dec)
            incl = numpy.copy(inclination).astype(numpy.float64)
            if opts.inclination_cosine_sampler:
              incl = numpy.arccos(incl)

            # use EXTREMELY many bits
            lnL = numpy.zeros(right_ascension.shape,dtype=numpy.float128)
            P.phi = right_ascension.astype(float)  # cast to float
            P.theta = dec #declination.astype(float)
            P.tref = float(fiducial_epoch)
            P.phiref = phi_orb.astype(float)
            P.incl = incl #inclination.astype(float)
            P.psi = psi.astype(float)
            P.dist = (distance* 1.e6 * lalsimutils.lsu_PC).astype(float) # luminosity distance

            lnL = factored_likelihood.DiscreteFactoredLogLikelihoodViaArrayVector(tvals,
                        P, lookupNKDict, rholmArrayDict, ctUArrayDict, ctVArrayDict,epochDict,Lmax=opts.l_max)
#            nEvals +=len(right_ascension)
            return numpy.exp(lnL-manual_avoid_overflow_logarithm)
        else:
            print( " Using CUDA GPU likelihood, if cupy available ")
            def likelihood_function(right_ascension, declination, phi_orb, inclination,
                psi, distance):
#                global nEvals
                tvals = xpy_default.linspace(-t_ref_wind,t_ref_wind,int((t_ref_wind)*2/P.deltaT))  # choose an array at the target sampling rate. P is inherited globally
                P.phi = xpy_default.asarray(right_ascension, dtype=np.float64)  # cast to float
                if opts.declination_cosine_sampler:
                  P.theta = numpy.pi/2 - xpy_default.arccos(xpy_default.asarray(declination,dtype=np.float64))
                else:
                  P.theta = xpy_default.asarray(declination,dtype=np.float64)
                P.tref = float(fiducial_epoch)
                P.phiref = xpy_default.asarray(phi_orb, dtype=np.float64)
                if opts.inclination_cosine_sampler:
                  P.incl = xpy_default.arccos(xpy_default.asarray(inclination, dtype=np.float64))
                else:
                  P.incl = xpy_default.asarray(inclination, dtype=np.float64)
                P.psi = xpy_default.asarray(psi, dtype=np.float64)
                P.dist = xpy_default.asarray(distance* 1.e6 * lalsimutils.lsu_PC,dtype=np.float64) # luminosity distance

                lnL = factored_likelihood.DiscreteFactoredLogLikelihoodViaArrayVectorNoLoop(tvals,
                        P, lookupNKDict, rholmArrayDict, ctUArrayDict, ctVArrayDict,epochDict,Lmax=opts.l_max,xpy=xpy_default)
#                nEvals +=len(right_ascension)
                return identity_convert(xpy_default.exp(lnL-manual_avoid_overflow_logarithm))

     else: # integrate over intrinsic variables. Right now those variables ahave HARDCODED NAMES, alas
        def likelihood_function(right_ascension, declination, phi_orb, inclination,
                psi, distance,q):
            dec = numpy.copy(declination).astype(numpy.float64)
            if opts.declination_cosine_sampler:
                dec = numpy.pi/2 - numpy.arccos(dec)
            incl = numpy.copy(inclination).astype(numpy.float64)
            if opts.inclination_cosine_sampler:
                incl = numpy.arccos(incl)

            lnL = numpy.zeros(len(right_ascension),dtype=numpy.float128)
#            i = 0
            tvals = numpy.linspace(-t_ref_wind,t_ref_wind,int((t_ref_wind)*2/P.deltaT))  # choose an array at the target sampling rate. P is inherited globally


#            t_start =lal.GPSTimeNow() 
            for ph, th, phr, ic, ps, di,qi in zip(right_ascension, dec,
                    phi_orb, incl, psi, distance,q):
                 # Reconstruct U,V using ROM fits.  PROBABLY should do this once for every q, rather than deep on the loop
                P.assign_param('q',qi)  # mass ratio
                rholms_intp_A, cross_terms_A, cross_terms_V_A, rholms_A, rest_A = factored_likelihood.ReconstructPrecomputedLikelihoodTermsROM(P, rest, rholms_intp, cross_terms, cross_terms_V, rholms,verbose=False)
                # proceed for rest
                P.phi = ph # right ascension
                P.theta = th # declination
                P.tref = fiducial_epoch  # see 'tvals', above
                P.phiref = phr # ref. orbital phase
                P.incl = ic # inclination
                P.psi = ps # polarization angle
                P.dist = di* 1.e6 * lalsimutils.lsu_PC # luminosity distance


                lnL[i] = factored_likelihood.FactoredLogLikelihoodTimeMarginalized(tvals,
                        P, rholms_intp_A, rholms_A, cross_terms_A, cross_terms_V_A,                   
                        opts.l_max,interpolate=opts.interpolate_time)
                if numpy.isnan(lnL[i]) or lnL[i]<-200:
                    lnL[i] = -200   # regularize  : a hack, for now, to deal with rare ROM problems. Only on the ROM logic fork
                i+=1
#            t_end =lal.GPSTimeNow() 
#            print " Cost per evaluation ", (t_end - t_start)/len(q)
#            print " Max lnL for this iteration ", numpy.max(lnL)
            return numpy.exp(lnL - manual_avoid_overflow_logarithm)

    if opts.sampler_method == "adaptive_cartesian_gpu":
      # reset sampling parameter as needed (distance, inclination but not sky location)
      # distance (and inclination) can conceivably correlate with mass, and we don't want to truncate in distance/inclination prematurely from history effects
      sampler.reset_sampling('distance')
      sampler.reset_sampling('inclination')

    # Integrate
    args = likelihood_function.__code__.co_varnames[:likelihood_function.__code__.co_argcount]
    print( " --------> Arguments ", args)
    res, var, neff, dict_return = sampler.integrate(likelihood_function, *unpinned_params, **pinned_params)


    # Report results
    if opts.output_file:
        fname_output_txt = opts.output_file +"_"+str(indx_event)+"_" + ".dat"
        m1 =P.m1/lal.MSUN_SI
        m2 =P.m2/lal.MSUN_SI
        if opts.sim_xml: 
            event_id = opts.event
        else:
            event_id = -1
        if opts.event == None:
            event_id = -1
        if opts.save_eccentricity:
            # output format when eccentricity is being used
            numpy.savetxt(fname_output_txt, numpy.array([[event_id, m1, m2, P.s1x, P.s1y, P.s1z, P.s2x, P.s2y, P.s2z,  P.eccentricity, numpy.log(res)+manual_avoid_overflow_logarithm, numpy.sqrt(var)/res,sampler.ntotal, neff ]]))  #dict_return["convergence_test_results"]["normal_integral]"
        elif not (P.lambda1>0 or P.lambda2>0):
            # output format when lambda is NOT used
            numpy.savetxt(fname_output_txt, numpy.array([[event_id, m1, m2, P.s1x, P.s1y, P.s1z, P.s2x, P.s2y, P.s2z,  numpy.log(res)+manual_avoid_overflow_logarithm, numpy.sqrt(var)/res,sampler.ntotal, neff ]]))  #dict_return["convergence_test_results"]["normal_integral]"
        else:
            # Alternative output format if lambda is active
            numpy.savetxt(fname_output_txt, numpy.array([[event_id, m1, m2, P.s1x, P.s1y, P.s1z, P.s2x, P.s2y, P.s2z,  P.lambda1, P.lambda2, numpy.log(res)+manual_avoid_overflow_logarithm, numpy.sqrt(var)/res,sampler.ntotal, neff ]]))  #dict_return["convergence_test_results"]["normal_integral]"

    # Comprehensive output (not yet provided)
    # Convert declination, inclination  parameters in sampler if needed
    if opts.save_samples and opts.output_file:
      import copy
      samples = copy.deepcopy(sampler._rvs)  # deep copy: avoid modifying structures and  having side effect on integrator, which loops over keys Expensive!
      if opts.inclination_cosine_sampler:
        samples["inclination"] = numpy.arccos(samples["inclination"].astype(numpy.float64))
      if opts.declination_cosine_sampler:
        samples["declination"] = numpy.pi/2 - numpy.arccos(samples["declination"].astype(numpy.float64))
      xmldoc = ligolw.Document()
      xmldoc.appendChild(ligolw.LIGO_LW())
      process.register_to_xmldoc(xmldoc, sys.argv[0], opts.__dict__)
      if not opts.time_marginalization:
            samples["t_ref"] += float(fiducial_epoch)
      else:
            samples["t_ref"] = float(fiducial_epoch)*numpy.ones(len(samples["psi"]))
      samples["polarization"] = samples["psi"]
      samples["coa_phase"] = samples["phi_orb"]
      if ("declination", "right_ascension") in sampler.params:
            samples["latitude"], samples["longitude"] = samples[("declination", "right_ascension")]
      else:
            samples["latitude"] = samples["declination"]
            samples["longitude"] = samples["right_ascension"]
      samples["loglikelihood"] = numpy.log(samples["integrand"]) + manual_avoid_overflow_logarithm  # export with consistent offset
      if not opts.rom_integrate_intrinsic:
            # ILE mode: insert fixed model parameters
            samples["mass1"] = numpy.ones(samples["psi"].shape)*m1 # opts.mass1
            samples["mass2"] = numpy.ones(samples["psi"].shape)*m2 # opts.mass2
            samples["spin1x"] =numpy.ones(samples["psi"].shape)*P.s1x
            samples["spin1y"] =numpy.ones(samples["psi"].shape)*P.s1y
            samples["spin1z"] =numpy.ones(samples["psi"].shape)*P.s1z
            samples["spin2x"] =numpy.ones(samples["psi"].shape)*P.s2x
            samples["spin2y"] =numpy.ones(samples["psi"].shape)*P.s2y
            samples["spin2z"] =numpy.ones(samples["psi"].shape)*P.s2z
      xmlutils.append_samples_to_xmldoc(xmldoc, samples)
        # Extra metadata
      dict_out={"mass1": m1, "mass2": m2, "alpha5":P.lambda1, "alpha6":P.lambda2, "event_duration": numpy.sqrt(var)/res, "ttotal": sampler.ntotal}
      if 'distance' in pinned_params:
            dict_out['distance'] = pinned_params["distance"]
      xmlutils.append_likelihood_result_to_xmldoc(xmldoc, numpy.log(res)+manual_avoid_overflow_logarithm, neff=neff, converged=dict_return["convergence_test_results"]["normal_integral"], **dict_out)
      fname_output_xml = opts.output_file +"_"+str(indx_event)+"_" + ".xml.gz"
      utils.write_filename(xmldoc, fname_output_xml, gz=True)


    # Clear sampler _rvs, to avoid side effects when called again
    for key in sampler._rvs.keys():
      sampler._rvs[key] = []

    return res


lnL_sofar = -np.inf
no_adapt_sky = False
for indx in numpy.arange(len(P_list)):
 try:
  res = analyze_event(P_list, indx, data_dict, psd_dict, fmax, opts)
  lnL_sofar = np.max([lnL_sofar,res])
  if opts.no_adapt_after_first and (not no_adapt_sky):
    if lnL_sofar > 20:  # Use absolute threshold.  Expect this will give modest sky localization.
      # remove right_ascension, declination from adaptive parameters
      params_adapt = sampler.adaptive
      params_adapt  = list(set(params_adapt) - set(['right_ascension','declination']))
      sampler.adaptive = params_adapt
      # Disable saving of the integrand -- saves on memory and will reduce speed. 
      # Note this needs to be done in a few places
      pinned_params.update({"force_no_adapt":True,"save_intg":False, "igrand_threshold_deltalnL":20})  # massively reduce memory usage in logic branch, don't save all. Highly redundant sequence to self-document all related logic branches
 except Exception as exception_failure:
  print( "  ===> FAILED ANALYSIS <==== ")
  print( exception_failure)
  if opts.internal_make_empty_file_on_error:
    fname_output_txt = opts.output_file +"_"+str(indx)+"_" + ".dat"
    open(fname_output_txt,'a').close()  # create empty file
  if opts.internal_hard_fail_on_error:
    sys.exit(1)
#  if len(exception_failure) >0:
#    if "CUBLAS" in exception_failure[0]:  # Hard fail if a cuda error!
#      sys.exit(1) 
#    if "CUDA" in exception_failure[0]: # Hard fail if a cuda error
#      sys.exit(1)
  print( " Probable reasons: SEOB nyquist or starting frequency limit or signal duration ")
  print( " Skipping the following binary! ")
  # Zero out extrinsic parameters -- these are CUDA-populated / meaningless, but could cause errors if populated
  P_list[indx].incl = P_list[indx].tref = P_list[indx].dist = P_list[indx].phiref = P_list[indx].psi =P_list[indx].theta = P_list[indx].phi =0
  P_list[indx].print_params()
  
