#!/usr/bin/env python

# Copyright (C) 2013-2019 Ian W. Harry, Alex Nitz, Marton Tapai,
#     Gareth Davies
#
# 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.
"""
Program for running multi-detector workflow analysis through coincidence and
then generate post-processing and plots.
"""
import pycbc
import pycbc.version
__author__  = "Alex Nitz <alex.nitz@ligo.org>"
__version__ = pycbc.version.git_verbose_msg
__date__    = pycbc.version.date
__program__ = "pycbc_offline"

import sys
import socket
import pycbc.events, pycbc.workflow as wf
import os, argparse, logging
import configparser as ConfigParser
from ligo import segments
import numpy, lal, datetime, itertools
from pycbc.results import create_versioning_page, static_table, layout
from pycbc.results.versioning import save_fig_with_metadata
from pycbc.results.metadata import html_escape


def symlink_path(f, path):
    if f is None:
        return
    try:
        os.symlink(f.storage_path, os.path.join(path, f.name))
    except OSError:
        pass

def symlink_result(f, rdir_path):
    symlink_path(f, rdir[rdir_path])


# Generator for producing ifo combinations
def ifo_combos(ifos):
    for i in range(2, len(ifos)+1):
        combinations = itertools.combinations(ifos, i)
        for ifocomb in combinations:
            yield ifocomb


# Log to the screen until we know where the output file is
logging.basicConfig(format='%(asctime)s:%(levelname)s : %(message)s',
    level=logging.INFO)

parser = argparse.ArgumentParser(description=__doc__[1:])
parser.add_argument('--version', action='version', version=__version__)
wf.add_workflow_command_line_group(parser)
wf.add_workflow_settings_cli(parser)
args = parser.parse_args()

# FIXME: opts.tags is currently unused here.

wf.makedir(args.output_dir)

container = wf.Workflow(args, name=args.workflow_name)
workflow = wf.Workflow(args, name=args.workflow_name + '-main')
finalize_workflow = wf.Workflow(args, name=args.workflow_name + '-finalization')
os.chdir(args.output_dir)

rdir = layout.SectionNumber('results', ['analysis_time',
                                 'detector_sensitivity',
                                 'data_quality',
                                 'single_triggers',
                                 'coincident_triggers',
                                 'injections',
                                 'search_sensitivity',
                                 'open_box_result',
                                 'workflow',
                                 ])

wf.makedir(rdir.base)
wf.makedir(rdir['workflow'])

wf_log_file = wf.File(workflow.ifos, 'workflow-log', workflow.analysis_time,
                      extension='.txt',
                      directory=rdir['workflow'])

logging.basicConfig(format='%(asctime)s:%(levelname)s : %(message)s',
                    filename=wf_log_file.storage_path,
                    level=logging.INFO,
                    filemode='w')

logfile = logging.FileHandler(filename=wf_log_file.storage_path,mode='w')
logfile.setLevel(logging.INFO)
formatter = logging.Formatter('%(asctime)s:%(levelname)s : %(message)s')
logfile.setFormatter(formatter)
logging.getLogger('').addHandler(logfile)
logging.info("Created log file %s" % wf_log_file.storage_path)

# put start / end time at top of summary page
time = workflow.analysis_time
s, e = int(time[0]), int(time[1])
s_utc = str(datetime.datetime(*lal.GPSToUTC(s)[0:6]))
e_utc = str(datetime.datetime(*lal.GPSToUTC(e)[0:6]))
time_str = '<center><p><b>GPS Interval [%s,%s). ' %(s,e)
time_str += 'UTC Interval %s - %s. ' %(s_utc, e_utc)
time_str += 'Interval duration = %.3f days.</b></p></center>'\
                                                         %(float(e-s)/86400.0,)
time_file = wf.File(workflow.ifos, 'time', workflow.analysis_time,
                                           extension='.html',
                                           directory=rdir.base)
kwds = { 'title' : 'Search Workflow Duration (Wall Clock Time)',
        'caption' : "Wall clock start and end times for this invocation of "
                    "the workflow. The command line button shows the "
                    "arguments used to invoke the workflow creation script.",
        'cmd' :' '.join(sys.argv), }
save_fig_with_metadata(time_str, time_file.storage_path, **kwds)

# Get segments and find the data locations
sci_seg_name = 'science'
science_seg_file = wf.get_segments_file(workflow, sci_seg_name, 'segments-science',
                                        rdir['analysis_time/segment_data'])

ssegs = {}
for ifo in workflow.ifos:
    ssegs[ifo] = science_seg_file.segment_dict["%s:science" % ifo]

hoft_tags=[]
if 'hoft' in workflow.cp.get_subsections('workflow-datafind'):
    hoft_tags=['hoft']
datafind_files, analyzable_file, analyzable_segs, analyzable_name = \
                                           wf.setup_datafind_workflow(workflow,
                                     ssegs, "datafind",
                                     seg_file=science_seg_file, tags=hoft_tags)

final_veto_name = 'vetoes'
final_veto_file = wf.get_segments_file(workflow, final_veto_name,
                                       'segments-vetoes',
                                       rdir['analysis_time/segment_data'])

# Get dq segments from veto definer and calculate data quality timeseries
dq_flag_name = 'dq_flag'
dq_segment_file = wf.get_flag_segments_file(workflow, dq_flag_name,
                                      'segments-dq',
                                      rdir['analysis_time/segment_data'])

# Template bank stuff
hdfbank = wf.setup_tmpltbank_workflow(workflow, analyzable_segs,
                                      datafind_files, output_dir="bank",
                                      return_format='hdf')
assert( len(hdfbank) == 1 )
hdfbank = hdfbank[0]

splitbank_files_fd = wf.setup_splittable_workflow(workflow, [hdfbank],
                                                  out_dir="bank",
                                                  tags=['full_data'])

bank_tags = []
if 'mass1_mass2' in workflow.cp.get_subsections('plot_bank'):
    bank_tags=['mass1_mass2']
bank_plot = wf.make_template_plot(workflow, hdfbank,
                                  rdir['coincident_triggers'],
                                  tags=bank_tags)


######################## Setup the FULL DATA run ##############################
output_dir = "full_data"

# setup the matchedfilter jobs
ind_insps = insps = wf.setup_matchedfltr_workflow(workflow, analyzable_segs,
                                   datafind_files, splitbank_files_fd,
                                   output_dir, tags=['full_data'])

insps = wf.merge_single_detector_hdf_files(workflow, hdfbank,
                                           insps, output_dir,
                                           tags=['full_data'])

# setup sngl trigger distribution fitting jobs
# 'statfiles' is list of files used in calculating coinc statistic
# 'dqfiles' is the subset of files containing data quality information
statfiles = []
dqfiles = []
dqfile_labels = []
dq_labels = workflow.cp.get_subsections('workflow-data_quality')
for dq_label in dq_labels:
    dq_label_files = wf.setup_dq_reranking(workflow, dq_label, insps, hdfbank,
                                       analyzable_segs, analyzable_file,
                                       dq_segment_file,
                                       output_dir='dq',
                                       tags=['full_data'])
    statfiles += dq_label_files
    dqfiles += dq_label_files
    dqfile_labels += len(dq_label_files) * [dq_label]
statfiles += wf.setup_trigger_fitting(workflow, insps, hdfbank,
                                      final_veto_file, final_veto_name,
                                      output_dir=output_dir,
                                      tags=['full_data'])

# Set up the multi-ifo coinc jobs
# final_bg_files contains coinc results using vetoes final_veto_files
# ifo_ids will store an (integer) index for each ifo in the precedence list
full_insps = insps
no_fg_exc_files = {}
ifo_ids = {}

# Get the ifo precedence values
ifo_precedence_list = workflow.cp.get_opt_tags('workflow-coincidence', 'timeslide-precedence', ['full_data'])
for ifo, _ in zip(*insps.categorize_by_attr('ifo')):
    ifo_ids[ifo] = ifo_precedence_list.index(ifo)

# Generate the possible detector combinations from 2 detectors
# up to the number of trigger files

if workflow.cp.has_option_tags('workflow-data_quality', 'no-coinc-veto',
                                tags=None):
    logging.info("no-coinc-veto option enabled, " +
                 "no longer passing veto segments to coinc jobs.")
    coinc_veto_file = None
else:
    coinc_veto_file = final_veto_file

for ifocomb in ifo_combos(ifo_ids.keys()):
    inspcomb = wf.select_files_by_ifo_combination(ifocomb, insps)
    pivot_ifo, fixed_ifo, ordered_ifo_list = wf.get_ordered_ifo_list(ifocomb,
                                                                     ifo_ids)

    # Create coinc tag, and set up the coinc job for the combination
    coinctag = '{}det'.format(len(ifocomb))
    ctagcomb = ['full_data', coinctag]

    bg_file = wf.setup_interval_coinc(
        workflow, hdfbank, inspcomb, statfiles, coinc_veto_file,
        final_veto_name, output_dir, pivot_ifo, fixed_ifo, tags=ctagcomb)

    # Optionally perform follow-up on triggers and rerank the candidates
    # Returns the input file if not enabled.
    no_fg_exc_files[ordered_ifo_list] = wf.rerank_coinc_followup(
                   workflow, bg_file, hdfbank, output_dir,
                   tags=ctagcomb)


if len(insps) == 2:
    final_bg_files = no_fg_exc_files
else:
    final_bg_files = {}
    for ifocomb in ifo_combos(ifo_ids.keys()):
        _, _, ordered_ifo_list = wf.get_ordered_ifo_list(ifocomb, ifo_ids)
        # Create coinc tag
        coinctag = '{}det'.format(len(ifocomb))
        ctagcomb = ['full_data', coinctag]
        other_ifo_keys = list(no_fg_exc_files.keys())
        other_ifo_keys.remove(ordered_ifo_list)
        other_bg_files = {ctype: no_fg_exc_files[ctype]
                          for ctype in other_ifo_keys}
        final_bg_files[ordered_ifo_list] = wf.setup_exclude_zerolag(
            workflow,
            no_fg_exc_files[ordered_ifo_list],
            wf.FileList(other_bg_files.values()),
            output_dir, ordered_ifo_list,
            tags=ctagcomb
        )

combined_bg_file = wf.setup_combine_statmap(
                                workflow,
                                wf.FileList(final_bg_files.values()),
                                wf.FileList([]),
                                output_dir,
                                tags=['full_data'])

censored_veto_name = 'closed_box'
censored_veto = wf.make_foreground_censored_veto(workflow,
                       combined_bg_file, final_veto_file, final_veto_name,
                       censored_veto_name, 'segments')

# Calculate the inspiral psds
psd_files = []
trig_generated_name = 'TRIGGERS_GENERATED'
trig_generated_segs = {}
data_analysed_name = 'DATA_ANALYSED'
data_analysed_segs = {}
insp_files_seg_dict = segments.segmentlistdict()

for ifo, files in zip(*ind_insps.categorize_by_attr('ifo')):
    trig_generated_segs[ifo] = segments.segmentlist([f.segment for f in files])
    data_analysed_segs[ifo] = \
        segments.segmentlist([f.metadata['data_seg'] for f in files])

    # Remove duplicates from splitbank
    trig_generated_segs[ifo] = \
        segments.segmentlist(set(trig_generated_segs[ifo]))
    data_analysed_segs[ifo] = \
        segments.segmentlist(set(data_analysed_segs[ifo]))

    insp_files_seg_dict[ifo + ":" + trig_generated_name] = \
                                                       trig_generated_segs[ifo]
    insp_files_seg_dict[ifo + ":" + data_analysed_name] = \
                                                        data_analysed_segs[ifo]

    if datafind_files:
        frame_files = datafind_files.find_output_with_ifo(ifo)
    else:
        frame_files = None
    psd_files += [wf.setup_psd_calculate(workflow, frame_files, ifo,
              data_analysed_segs[ifo], data_analysed_name, 'psds')]

insp_files_seg_file = wf.SegFile.from_segment_list_dict('INSP_SEGMENTS',
                 insp_files_seg_dict, valid_segment=workflow.analysis_time,
                 extension='xml', directory=rdir['analysis_time/segment_data'])

################### Range, spectrum and segments plots #######################

s = wf.make_spectrum_plot(workflow, psd_files, rdir['detector_sensitivity'])
r = wf.make_range_plot(workflow, psd_files, rdir['detector_sensitivity'],
                       require='summ')
r2 = wf.make_range_plot(workflow, psd_files, rdir['detector_sensitivity'],
                        exclude='summ')

det_summ = [(s, r[0] if len(r) != 0 else None)]
layout.two_column_layout(rdir['detector_sensitivity'],
                         det_summ + list(layout.grouper(r2, 2)))

# do plotting of segments / veto times
wf.make_segments_plot(workflow, [insp_files_seg_file],
                      rdir['analysis_time/segments'],
                      tags=[trig_generated_name])
wf.make_gating_plot(workflow, full_insps, rdir['analysis_time/gating'],
                    tags=['full_data'])

# make segment table and plot for summary page
curr_files = [science_seg_file, analyzable_file,
              insp_files_seg_file]
curr_names = [sci_seg_name, analyzable_name,
              trig_generated_name]
seg_summ_table = wf.make_seg_table\
    (workflow, curr_files, curr_names, rdir['analysis_time/segments'],
     ['SUMMARY'], title_text='Input and output',
     description='This shows the total amount of input data, analyzable data, '
                 'and the time for which triggers are produced.')
seg_summ_plot = wf.make_seg_plot(workflow, curr_files,
                                rdir['analysis_time/segments'],
                                curr_names, ['SUMMARY'])

curr_files = [insp_files_seg_file] + [final_veto_file]
# Add in singular veto files
curr_files = curr_files + [science_seg_file]
curr_names = [trig_generated_name + '&' + final_veto_name]
# And SCIENCE - CAT 1 vetoes explicitly.
curr_names += [sci_seg_name + '&' + 'VETO_CAT1']

veto_summ_table = wf.make_seg_table\
    (workflow, curr_files, curr_names, rdir['analysis_time/segments'],
     ['VETO_SUMMARY'], title_text='Time removed by vetoes',
     description='This shows the time removed from the output time by the '
                 'vetoes applied to the triggers.')

# make veto definer table
vetodef_table = wf.make_veto_table(workflow, rdir['analysis_time/veto_definer'])
if vetodef_table is not None:
    layout.single_layout(rdir['analysis_time/veto_definer'], ([vetodef_table]))

#################### Plotting on FULL_DATA results ##########################
##################### SINGLES plots first ###################################

snrchi = wf.make_snrchi_plot(workflow, insps, censored_veto,
                            'closed_box', rdir['single_triggers'],
                            tags=['full_data'])
layout.group_layout(rdir['single_triggers'], snrchi)

hist_summ = []
for insp in full_insps:
    outdir = rdir['single_triggers/%s_binned_triggers' % insp.ifo]
    wf.make_singles_plot(workflow, [insp], hdfbank, censored_veto,
                         'closed_box', outdir, tags=['full_data'])
    # make non-summary hists using the bank file
    # currently, none of these are made
    outdir = rdir['single_triggers/%s_trigger_histograms' % insp.ifo]
    wf.make_single_hist(workflow, insp, censored_veto, 'closed_box', outdir,
                        bank_file=hdfbank, exclude='summ',
                        tags=['full_data'])
    # make summary hists for all templates together
    # currently, 2 per ifo: snr and newsnr
    outdir = rdir['single_triggers/%s_trigger_histograms' % insp.ifo]
    allhists = wf.make_single_hist(workflow, insp, censored_veto, 'closed_box',
                                   outdir, require='summ', tags=['full_data'])
    # make hists of newsnr split up by parameter
    # currently, 1 per ifo split by template duration
    outdir = rdir['single_triggers/%s_binned_histograms' % insp.ifo]
    binhists = wf.make_binned_hist(workflow, insp, censored_veto,
                                   'closed_box', outdir, hdfbank,
                                   tags=['full_data'])
    # put raw SNR and binned newsnr hist in summary
    hist_summ += list(layout.grouper([allhists[0], binhists[0]], 2))

if workflow.cp.has_option_tags('workflow-matchedfilter',
                               'plot-throughput', tags=['full_data']):
    wf.make_throughput_plot(workflow, full_insps, rdir['workflow/throughput'],
                            tags=['full_data'])

# Run minifollowups on loudest sngl detector events
excl_subsecs = set([])

for insp_file in full_insps:
    for tag in insp_file.tags:
        excl_subsecs.add(tag)

for insp_file in full_insps:
    curr_ifo = insp_file.ifo
    for subsec in workflow.cp.get_subsections('workflow-sngl_minifollowups'):
        if subsec in excl_subsecs:
            continue
        sec_name = 'workflow-sngl_minifollowups-{}'.format(subsec)
        dir_str = workflow.cp.get(sec_name, 'section-header')
        currdir = rdir['single_triggers/{}_{}'.format(curr_ifo, dir_str)]
        wf.setup_single_det_minifollowups\
            (workflow, insp_file, hdfbank, insp_files_seg_file,
             data_analysed_name, trig_generated_name, 'daxes', currdir,
             veto_file=censored_veto, veto_segment_name='closed_box',
             tags=insp_file.tags + [subsec])

##################### COINC FULL_DATA plots ###################################

# Main results with combined file (we mix open and closed box here, but
# separate them in the result page)

# FIXME: COMMENTED OUT JOBS ARE FAILING ... NEED FIXING!!
#        (Currently that's most of the jobs :-( )
#snrifar = wf.make_snrifar_plot(workflow, combined_bg_file,
#                               rdir['open_box_result'],
#                               tags=combined_bg_file.tags)
#snrifar_cb = wf.make_snrifar_plot(workflow, combined_bg_file,
#                                  rdir['coincident_triggers'], closed_box=True,
#                                  tags=combined_bg_file.tags + ['closed'])
#ratehist = wf.make_snrratehist_plot(workflow, combined_bg_file,
#                                    rdir['open_box_result'],
#                                    tags=combined_bg_file.tags)
#snrifar_ifar = wf.make_snrifar_plot(workflow, combined_bg_file,
#                                    rdir['open_box_result/significance'],
#                                    cumulative=False,
#                                    tags=combined_bg_file.tags + ['ifar'])
ifar_ob = wf.make_ifar_plot(workflow, combined_bg_file,
                                    rdir['open_box_result'],
                                    tags=combined_bg_file.tags + ['open_box'],
                                    executable='page_ifar_catalog')
#ifar_cb = wf.make_ifar_plot(workflow, combined_bg_file,
#                                    rdir['coincident_triggers'],
#                                    tags=combined_bg_file.tags + ['closed_box'],
#                                    executable='page_ifar_catalog')
table = wf.make_foreground_table(workflow, combined_bg_file,
                                 hdfbank, rdir['open_box_result'],
                                 singles=insps, extension='.html',
                                 tags=combined_bg_file.tags)

fore_xmlall = wf.make_foreground_table(workflow, combined_bg_file,
                    hdfbank, rdir['open_box_result'], singles=insps,
                    extension='.xml', tags=["xmlall"])
fore_xmlloudest = wf.make_foreground_table(workflow, combined_bg_file,
                    hdfbank, rdir['open_box_result'], singles=insps,
                    extension='.xml', tags=["xmlloudest"])

#symlink_result(snrifar, 'open_box_result/significance')
#symlink_result(ratehist, 'open_box_result/significance')
symlink_result(table, 'open_box_result/significance')

# Set html pages
main_page = [(ifar_ob,), (table, )]
layout.two_column_layout(rdir['open_box_result'], main_page)

#detailed_page = [(snrifar, ratehist), (snrifar_ifar, ifar_ob), (table,)]
#layout.two_column_layout(rdir['open_box_result/significance'], detailed_page)

closed_page = [(bank_plot,)]
layout.two_column_layout(rdir['coincident_triggers'], closed_page)

# run minifollowups on the output of the loudest events
mfup_dir_fg = rdir['open_box_result/loudest_events_followup']
mfup_dir_bg = rdir['coincident_triggers/loudest_background_followup']
wf.setup_foreground_minifollowups(workflow, combined_bg_file,
                                  full_insps, hdfbank, insp_files_seg_file,
                                  data_analysed_name, trig_generated_name,
                                  'daxes', mfup_dir_fg,
                                  tags=combined_bg_file.tags + ['foreground'])

wf.setup_foreground_minifollowups(workflow, combined_bg_file,
                                  full_insps, hdfbank, insp_files_seg_file,
                                  data_analysed_name, trig_generated_name,
                                  'daxes', mfup_dir_bg,
                                  tags=combined_bg_file.tags + ['background'])

# Sub-pages for each ifo combination

snrifar_summ = []
for key in final_bg_files:
    # FIXME: Stop obfuscating this file!
    bg_file = final_bg_files[key]
    open_dir = rdir['open_box_result/{}_coincidences'.format(key)]
    closed_dir = rdir['coincident_triggers/{}_coincidences'.format(key)]
    snrifar = wf.make_snrifar_plot(workflow, bg_file, open_dir,
                                   tags=bg_file.tags)
    snrifar_cb = wf.make_snrifar_plot(workflow, bg_file, closed_dir,
                                      closed_box=True,
                                      tags=bg_file.tags + ['closed'])
    ratehist = wf.make_snrratehist_plot(workflow, bg_file, open_dir,
                                        tags=bg_file.tags)
    snrifar_ifar = wf.make_snrifar_plot(workflow, bg_file, open_dir,
                                        cumulative=False,
                                        tags=bg_file.tags + ['ifar'])
    ifar_ob = wf.make_ifar_plot(workflow, bg_file, open_dir,
                                tags=bg_file.tags + ['open_box'])
    ifar_cb = wf.make_ifar_plot(workflow, bg_file, closed_dir,
                                tags=bg_file.tags + ['closed_box'])
    table = wf.make_foreground_table(workflow, bg_file, hdfbank, open_dir,
                                     singles=insps, extension='.html',
                                     tags=bg_file.tags)

    detailed_page = [(snrifar, ratehist), (snrifar_ifar, ifar_ob), (table,)]
    layout.two_column_layout(open_dir, detailed_page)

    closed_page = [(snrifar_cb, ifar_cb)]
    snrifar_summ += closed_page
    layout.two_column_layout(closed_dir, closed_page)

#################### Plotting of data quality results #########################

# DQ log likelihood plots
for dqf,dql in zip(dqfiles,dqfile_labels):
    outdir = rdir[f'data_quality/{dqf.ifo}_trigger_rates_{dql}']
    wf.make_dq_trigger_rate_plot(workflow, [dqf], outdir, tags=[dql,'full_data'])

    outdir = rdir[f'data_quality/{dqf.ifo}_percentiles_{dql}']
    wf.make_dq_percentile_plot(workflow, [dqf], outdir, tags=[dql,'full_data'])

# DQ background binning plot
if workflow.cp.has_option_tags('workflow-data_quality', 'background-bins',
                                tags=None):
    background_bins = workflow.cp.get_opt_tags('bin_trigger_rates_dq',
                                               'background-bins',tags=None)
    bank_tags = workflow.cp.get_subsections('plot_bank')
    if len(bank_tags) == 0: bank_tags=['bank']
    for b_tag in bank_tags:
        wf.make_template_plot(workflow, hdfbank, rdir['data_quality'],
                              bins=background_bins, tags=['dq_bins',b_tag])

############################## Setup the injection runs #######################

splitbank_files_inj = wf.setup_splittable_workflow(workflow, [hdfbank],
                                                   out_dir="bank",
                                                   tags=['injections'])

# setup the injection files
inj_files_base, inj_tags = wf.setup_injection_workflow(workflow,
                                                  output_dir="inj_files")
inj_files = []
for inj_file, tag in zip(inj_files_base, inj_tags):
    inj_files.append(wf.inj_to_hdf(workflow, inj_file, 'inj_files', [tag]))

inj_coincs = wf.FileList()

found_inj_dict ={}
insps_dict = {}

files_for_combined_injfind = []
for inj_file, tag in zip(inj_files, inj_tags):
    ctags = [tag, 'injections']
    output_dir = '%s_coinc' % tag

    if workflow.cp.has_option_tags('workflow-injections',
                                   'compute-optimal-snr', tags=ctags):
        optimal_snr_file = wf.compute_inj_optimal_snr(
                workflow, inj_file, psd_files, 'inj_files', tags=ctags)
        file_for_injfind = optimal_snr_file
    else:
        file_for_injfind = inj_file

    files_for_combined_injfind.append(file_for_injfind)

    # setup the matchedfilter jobs
    insps = wf.setup_matchedfltr_workflow(workflow, analyzable_segs,
                                         datafind_files, splitbank_files_inj,
                                         output_dir, tags=ctags,
                                         injection_file=inj_file)

    insps = wf.merge_single_detector_hdf_files(workflow, hdfbank,
                                               insps, output_dir, tags=ctags)
    # coincidence for injections
    inj_coinc = {}
    for ifocomb in ifo_combos(ifo_ids.keys()):
        inspcomb = wf.select_files_by_ifo_combination(ifocomb, insps)
        pivot_ifo, fixed_ifo, ordered_ifo_list = \
            wf.get_ordered_ifo_list(ifocomb, ifo_ids)

        # Create coinc tag, and set up the coinc job for the combination
        coinctag = '{}det'.format(len(ifocomb))
        ctagcomb = [tag, 'injections', coinctag]
        curr_out = wf.setup_interval_coinc_inj(
            workflow,
            hdfbank,
            full_insps,
            inspcomb,
            statfiles,
            final_bg_files[ordered_ifo_list],
            coinc_veto_file,
            final_veto_name,
            output_dir,
            pivot_ifo,
            fixed_ifo,
            tags=ctagcomb
        )

        # Rerank events
        curr_out = wf.rerank_coinc_followup(
            workflow,
            curr_out,
            hdfbank,
            output_dir,
            tags=ctagcomb,
            injection_file=inj_file,
            ranking_file=final_bg_files[ordered_ifo_list]
        )

        inj_coinc[ordered_ifo_list] = curr_out

    combctags = [tag, 'injections']
    final_inj_bg_file_list = wf.FileList(inj_coinc.values())
    combined_inj_bg_file = wf.setup_combine_statmap(
        workflow,
        final_inj_bg_file_list,
        wf.FileList(final_bg_files.values()),
        output_dir,
        tags=combctags
    )

    found_inj = wf.find_injections_in_hdf_coinc(
        workflow,
        [combined_inj_bg_file],
        [file_for_injfind],
        censored_veto,
        censored_veto_name,
        output_dir,
        tags=combctags)

    inj_coincs += [combined_inj_bg_file]

    # Set files for plots
    found_inj_dict[tag] = found_inj
    insps_dict[tag] = insps

# And the combined INJFIND file
if len(files_for_combined_injfind) > 0:
    found_inj_comb = wf.find_injections_in_hdf_coinc(
        workflow,
        inj_coincs,
        files_for_combined_injfind,
        censored_veto,
        censored_veto_name,
        'allinj',
        tags=['all','injections']
    )

############################ Injection plots #################################
# Per injection run plots
for inj_file, tag in zip(inj_files, inj_tags):
    found_inj = found_inj_dict[tag]
    insps = insps_dict[tag]
    injdir = rdir['injections/%s' % tag]
    sensdir = rdir['search_sensitivity/%s' % tag]

    #foundmissed/sensitivity plots
    s = wf.make_sensitivity_plot(workflow, found_inj, sensdir,
                                 exclude=['all', 'summ'], require='sub',
                                 tags=[tag])
    f = wf.make_foundmissed_plot(workflow, found_inj, injdir,
                                 exclude=['all', 'summ'], require='sub',
                                 tags=[tag])

    found_table = wf.make_inj_table(workflow, found_inj, injdir,
                                    singles=insps, tags=[tag, 'found'])
    missed_table = wf.make_inj_table(workflow, found_inj, injdir,
                                     missed=True, tags=[tag, 'missed'])

    for ifocomb in ifo_combos(ifo_ids.keys()):
        inspcomb = wf.select_files_by_ifo_combination(ifocomb, insps)
        fdinspcmb = wf.select_files_by_ifo_combination(ifocomb, full_insps)
        _, _, ordered_ifo_list = wf.get_ordered_ifo_list(ifocomb, ifo_ids)

        for inj_insp, trig_insp in zip(inspcomb, fdinspcmb):
            final_bg_file = final_bg_files[ordered_ifo_list]
            curr_tags = [tag, ordered_ifo_list]
            f += wf.make_coinc_snrchi_plot(workflow, found_inj, inj_insp,
                                           final_bg_file, trig_insp,
                                           injdir, tags=curr_tags)

    # Make pages from plots
    inj_layout = list(layout.grouper(f, 2)) + [(found_table,), (missed_table,)]
    layout.two_column_layout(rdir['injections/%s' % tag], inj_layout)

    if len(s) > 0:
        layout.group_layout(rdir['search_sensitivity/%s' % tag], s)

    # Minifollowups
    curr_dir_nam = 'injections/followup_of_{}'.format(tag)
    if workflow.cp.has_option_tags('workflow-injection_minifollowups',
                                   'subsection-suffix', tags=[tag]):
        suf_str = workflow.cp.get_opt_tags('workflow-injection_minifollowups',
                                           'subsection-suffix', tags=[tag])
        curr_dir_nam += '_' + suf_str
    currdir = rdir[curr_dir_nam]
    wf.setup_injection_minifollowups(workflow, found_inj, inj_file,
                                     insps, hdfbank, insp_files_seg_file,
                                     data_analysed_name, trig_generated_name,
                                     'daxes', currdir, tags=[tag])

    # If option given, make throughput plots
    if workflow.cp.has_option_tags('workflow-matchedfilter',
                                   'plot-throughput', tags=[tag]):
        wf.make_throughput_plot(workflow, insps, rdir['workflow/throughput'],
                                tags=[tag])


######################## Make combined injection plots ##########################
if len(files_for_combined_injfind) > 0:
    sen_all = wf.make_sensitivity_plot(workflow, found_inj_comb,
                                       rdir['search_sensitivity'],
                                       require='all')
    layout.group_layout(rdir['search_sensitivity'], sen_all)
    inj_all = wf.make_foundmissed_plot(workflow, found_inj_comb,
                                       rdir['injections'],
                                       require='all')
    layout.group_layout(rdir['injections'], inj_all)

    # Make summary page foundmissed and sensitivity plot
    sen_s = wf.make_sensitivity_plot(workflow, found_inj_comb,
                                     rdir['search_sensitivity'],
                                     require='summ')
    inj_s = wf.make_foundmissed_plot(workflow, found_inj_comb,
                                     rdir['injections'],
                                     require='summ')
    inj_summ = list(layout.grouper(inj_s + sen_s, 2))


# Make analysis time summary
analysis_time_summ = [time_file, seg_summ_plot]
for f in analysis_time_summ:
    symlink_result(f, 'analysis_time')
layout.single_layout(rdir['analysis_time'], (analysis_time_summ))


########################## Make full summary ####################################
if len(files_for_combined_injfind) > 0:
    summ = ([(time_file,)] + [(seg_summ_plot,)] +
            [(seg_summ_table, veto_summ_table)] + det_summ + hist_summ +
            [(bank_plot,)] + inj_summ + snrifar_summ)

else:
    summ = ([(time_file,)] + [(seg_summ_plot,)] +
            [(seg_summ_table, veto_summ_table)] + det_summ + hist_summ +
            [(bank_plot,)] + snrifar_summ)

for row in summ:
    for f in row:
        symlink_path(f, rdir.base)
layout.two_column_layout(rdir.base, summ)


# Save global config file to results directory
base = rdir['workflow/configuration']
wf.makedir(base)
ini_file_path = os.path.join(base, 'configuration.ini')
with open(ini_file_path, 'w') as ini_fh:
    container.cp.write(ini_fh)
ini_file = wf.FileList([wf.File(workflow.ifos, '', workflow.analysis_time,
                        file_url='file://' + ini_file_path)])
layout.single_layout(base, ini_file)


# Create versioning information
create_versioning_page(rdir['workflow/version'], container.cp)


############################ Finalization ####################################

# Create the final log file
log_file_html = wf.File(workflow.ifos, 'WORKFLOW-LOG', workflow.analysis_time,
                        extension='.html', directory=rdir['workflow'])

gen_file_html = wf.File(workflow.ifos, 'WORKFLOW-GEN', workflow.analysis_time,
                        extension='.html', directory=rdir['workflow'])

# Create a page to contain a dashboard link
dashboard_file = wf.File(workflow.ifos, 'DASHBOARD', workflow.analysis_time,
                         extension='.html', directory=rdir['workflow'])
dashboard_str = """<center><p style="font-size:20px"><b><a href="PEGASUS_DASHBOARD_URL" target="_blank">Pegasus Dashboard Page</a></b></p></center>"""
kwds = { 'title' : "Pegasus Dashboard",
         'caption' : "Link to Pegasus Dashboard",
         'cmd' : "PYCBC_SUBMIT_DAX_ARGV", }
save_fig_with_metadata(dashboard_str, dashboard_file.storage_path, **kwds)

# Create pages for the submission script to write data
wf.makedir(rdir['workflow/dax'])
wf.makedir(rdir['workflow/input_map'])
wf.makedir(rdir['workflow/output_map'])
wf.makedir(rdir['workflow/planning'])

wf.make_results_web_page(finalize_workflow, os.path.join(os.getcwd(),
                         rdir.base))

container += workflow
container += finalize_workflow

container.add_subworkflow_dependancy(workflow, finalize_workflow)

container.save()

# Protect the open box results folder
os.chmod(rdir['open_box_result'], 0o0700)

logging.info("Written dax.")

# Close the log and flush to the html file
logging.shutdown()
with open (wf_log_file.storage_path, "r") as logfile:
    logdata=logfile.read()
log_str = """
<p>Workflow generation script created workflow in output directory: %s</p>
<p>Workflow name is: %s</p>
<p>Workflow generation script run on host: %s</p>
<pre>%s</pre>
""" % (os.getcwd(), args.workflow_name, socket.gethostname(), logdata)
kwds = { 'title' : 'Workflow Generation Log',
         'caption' : "Log of the workflow script %s" % sys.argv[0],
         'cmd' :' '.join(sys.argv), }
save_fig_with_metadata(log_str, log_file_html.storage_path, **kwds)

# Add the command line used to a specific file
args_to_output = [sys.argv[0]]
for arg in sys.argv[1:]:
    if arg.startswith('--'):
        # This is an option, add tab
        args_to_output.append('  ' + arg)
    else:
        # This is a parameter, add two tabs
        args_to_output.append('    ' + arg)

gen_str = '<pre>' + ' \\\n'.join(args_to_output) + '</pre>'
kwds = { 'title' : 'Workflow Generation Command',
         'caption' : "Command used to generate the workflow.",
         'cmd' :' '.join(sys.argv), }
save_fig_with_metadata(gen_str, gen_file_html.storage_path, **kwds)
layout.single_layout(rdir['workflow'], ([dashboard_file, gen_file_html, log_file_html]))
