estrain@0: #!/usr/bin/env python3 estrain@0: estrain@0: import sys estrain@0: import time estrain@0: import random estrain@0: import os estrain@0: import subprocess estrain@0: import gzip estrain@0: import io estrain@0: import pickle estrain@0: import argparse estrain@0: import itertools estrain@0: import math estrain@0: from distutils.version import LooseVersion estrain@0: from distutils.spawn import find_executable estrain@0: estrain@0: try: estrain@0: from .version import SeqSero2_version estrain@0: except Exception: #ImportError estrain@0: from version import SeqSero2_version estrain@0: estrain@0: ### SeqSero Kmer estrain@0: def parse_args(): estrain@0: "Parse the input arguments, use '-h' for help." estrain@0: parser = argparse.ArgumentParser(usage='SeqSero2_package.py -t -m -i [-d ] [-p ] [-b ]\n\nDevelopper: Shaokang Zhang (zskzsk@uga.edu), Hendrik C Den-Bakker (Hendrik.DenBakker@uga.edu) and Xiangyu Deng (xdeng@uga.edu)\n\nContact email:seqsero@gmail.com\n\nVersion: v1.0.2')#add "-m " in future estrain@0: parser.add_argument("-i",nargs="+",help=": path/to/input_data",type=os.path.abspath) ### ed_SL_05282019: add 'type=os.path.abspath' to generate absolute path of input data. estrain@0: parser.add_argument("-t",choices=['1','2','3','4','5','6'],help=": '1' for interleaved paired-end reads, '2' for separated paired-end reads, '3' for single reads, '4' for genome assembly, '5' for nanopore fasta, '6' for nanopore fastq") estrain@0: parser.add_argument("-b",choices=['sam','mem'],default="mem",help=": algorithms for bwa mapping for allele mode; 'mem' for mem, 'sam' for samse/sampe; default=mem; optional; for now we only optimized for default 'mem' mode") estrain@0: parser.add_argument("-p",default="1",help=": number of threads for allele mode, if p >4, only 4 threads will be used for assembly since the amount of extracted reads is small, default=1") estrain@0: parser.add_argument("-m",choices=['k','a'],default="a",help=": which workflow to apply, 'a'(raw reads allele micro-assembly), 'k'(raw reads and genome assembly k-mer), default=a") estrain@0: parser.add_argument("-d",help=": output directory name, if not set, the output directory would be 'SeqSero_result_'+time stamp+one random number") estrain@0: parser.add_argument("-c",action="store_true",help=": if '-c' was flagged, SeqSero2 will only output serotype prediction without the directory containing log files") estrain@0: parser.add_argument("--check",action="store_true",help=": use '--check' flag to check the required dependencies") estrain@0: parser.add_argument('-v', '--version', action='version', version='%(prog)s ' + SeqSero2_version) estrain@0: return parser.parse_args() estrain@0: estrain@0: ### ed_SL_05282019: check paths of dependencies estrain@0: check_dependencies = parse_args().check estrain@0: dependencies = ['bwa','samtools','blastn','fastq-dump','spades.py','bedtools','SalmID.py'] estrain@0: if check_dependencies: estrain@0: for item in dependencies: estrain@0: ext_path = find_executable(item) estrain@0: if ext_path is not None: estrain@0: print ("Using "+item+" - "+ext_path) estrain@0: else: estrain@0: print ("ERROR: can not find "+item+" in PATH") estrain@0: sys.exit() estrain@0: ### end of --check estrain@0: estrain@0: def reverse_complement(sequence): estrain@0: complement = { estrain@0: 'A': 'T', estrain@0: 'C': 'G', estrain@0: 'G': 'C', estrain@0: 'T': 'A', estrain@0: 'N': 'N', estrain@0: 'M': 'K', estrain@0: 'R': 'Y', estrain@0: 'W': 'W', estrain@0: 'S': 'S', estrain@0: 'Y': 'R', estrain@0: 'K': 'M', estrain@0: 'V': 'B', estrain@0: 'H': 'D', estrain@0: 'D': 'H', estrain@0: 'B': 'V' estrain@0: } estrain@0: return "".join(complement[base] for base in reversed(sequence)) estrain@0: estrain@0: estrain@0: def createKmerDict_reads(list_of_strings, kmer): estrain@0: kmer_table = {} estrain@0: for string in list_of_strings: estrain@0: sequence = string.strip('\n') estrain@0: for i in range(len(sequence) - kmer + 1): estrain@0: new_mer = sequence[i:i + kmer].upper() estrain@0: new_mer_rc = reverse_complement(new_mer) estrain@0: if new_mer in kmer_table: estrain@0: kmer_table[new_mer.upper()] += 1 estrain@0: else: estrain@0: kmer_table[new_mer.upper()] = 1 estrain@0: if new_mer_rc in kmer_table: estrain@0: kmer_table[new_mer_rc.upper()] += 1 estrain@0: else: estrain@0: kmer_table[new_mer_rc.upper()] = 1 estrain@0: return kmer_table estrain@0: estrain@0: estrain@0: def multifasta_dict(multifasta): estrain@0: multifasta_list = [ estrain@0: line.strip() for line in open(multifasta, 'r') if len(line.strip()) > 0 estrain@0: ] estrain@0: headers = [i for i in multifasta_list if i[0] == '>'] estrain@0: multifasta_dict = {} estrain@0: for h in headers: estrain@0: start = multifasta_list.index(h) estrain@0: for element in multifasta_list[start + 1:]: estrain@0: if element[0] == '>': estrain@0: break estrain@0: else: estrain@0: if h[1:] in multifasta_dict: estrain@0: multifasta_dict[h[1:]] += element estrain@0: else: estrain@0: multifasta_dict[h[1:]] = element estrain@0: return multifasta_dict estrain@0: estrain@0: estrain@0: def multifasta_single_string(multifasta): estrain@0: multifasta_list = [ estrain@0: line.strip() for line in open(multifasta, 'r') estrain@0: if (len(line.strip()) > 0) and (line.strip()[0] != '>') estrain@0: ] estrain@0: return ''.join(multifasta_list) estrain@0: estrain@0: estrain@0: def chunk_a_long_sequence(long_sequence, chunk_size=60): estrain@0: chunk_list = [] estrain@0: steps = len(long_sequence) // 60 #how many chunks estrain@0: for i in range(steps): estrain@0: chunk_list.append(long_sequence[i * chunk_size:(i + 1) * chunk_size]) estrain@0: chunk_list.append(long_sequence[steps * chunk_size:len(long_sequence)]) estrain@0: return chunk_list estrain@0: estrain@0: estrain@0: def target_multifasta_kmerizer(multifasta, k, kmerDict): estrain@0: forward_length = 300 #if find the target, put forward 300 bases estrain@0: reverse_length = 2200 #if find the target, put backward 2200 bases estrain@0: chunk_size = 60 #it will firstly chunk the single long sequence to multiple smaller sequences, it controls the size of those smaller sequences estrain@0: target_mers = [] estrain@0: long_single_string = multifasta_single_string(multifasta) estrain@0: multifasta_list = chunk_a_long_sequence(long_single_string, chunk_size) estrain@0: unit_length = len(multifasta_list[0]) estrain@0: forward_lines = int(forward_length / unit_length) + 1 estrain@0: reverse_lines = int(forward_length / unit_length) + 1 estrain@0: start_num = 0 estrain@0: end_num = 0 estrain@0: for i in range(len(multifasta_list)): estrain@0: if i not in range(start_num, end_num): #avoid computational repetition estrain@0: line = multifasta_list[i] estrain@0: start = int((len(line) - k) // 2) estrain@0: s1 = line[start:k + start] estrain@0: if s1 in kmerDict: #detect it is a potential read or not (use the middle part) estrain@0: if i - forward_lines >= 0: estrain@0: start_num = i - forward_lines estrain@0: else: estrain@0: start_num = 0 estrain@0: if i + reverse_lines <= len(multifasta_list) - 1: estrain@0: end_num = i + reverse_lines estrain@0: else: estrain@0: end_num = len(multifasta_list) - 1 estrain@0: target_list = [ estrain@0: x.strip() for x in multifasta_list[start_num:end_num] estrain@0: ] estrain@0: target_line = "".join(target_list) estrain@0: target_mers += [ estrain@0: k1 for k1 in createKmerDict_reads([str(target_line)], k) estrain@0: ] ##changed k to k1, just want to avoid the mixes of this "k" (kmer) to the "k" above (kmer length) estrain@0: else: estrain@0: pass estrain@0: return set(target_mers) estrain@0: estrain@0: estrain@0: def target_read_kmerizer(file, k, kmerDict): estrain@0: i = 1 estrain@0: n_reads = 0 estrain@0: total_coverage = 0 estrain@0: target_mers = [] estrain@0: if file.endswith(".gz"): estrain@0: file_content = io.BufferedReader(gzip.open(file)) estrain@0: else: estrain@0: file_content = open(file, "r").readlines() estrain@0: for line in file_content: estrain@0: start = int((len(line) - k) // 2) estrain@0: if i % 4 == 2: estrain@0: if file.endswith(".gz"): estrain@0: s1 = line[start:k + start].decode() estrain@0: line = line.decode() estrain@0: else: estrain@0: s1 = line[start:k + start] estrain@0: if s1 in kmerDict: #detect it is a potential read or not (use the middle part) estrain@0: n_reads += 1 estrain@0: total_coverage += len(line) estrain@0: target_mers += [ estrain@0: k1 for k1 in createKmerDict_reads([str(line)], k) estrain@0: ] #changed k to k1, just want to avoid the mixes of this "k" (kmer) to the "k" above (kmer length) estrain@0: i += 1 estrain@0: if total_coverage >= 4000000: estrain@0: break estrain@0: return set(target_mers) estrain@0: estrain@0: estrain@0: def minion_fasta_kmerizer(file, k, kmerDict): estrain@0: i = 1 estrain@0: n_reads = 0 estrain@0: total_coverage = 0 estrain@0: target_mers = {} estrain@0: for line in open(file): estrain@0: if i % 2 == 0: estrain@0: for kmer, rc_kmer in kmers(line.strip().upper(), k): estrain@0: if (kmer in kmerDict) or (rc_kmer in kmerDict): estrain@0: if kmer in target_mers: estrain@0: target_mers[kmer] += 1 estrain@0: else: estrain@0: target_mers[kmer] = 1 estrain@0: if rc_kmer in target_mers: estrain@0: target_mers[rc_kmer] += 1 estrain@0: else: estrain@0: target_mers[rc_kmer] = 1 estrain@0: i += 1 estrain@0: return set([h for h in target_mers]) estrain@0: estrain@0: estrain@0: def minion_fastq_kmerizer(file, k, kmerDict): estrain@0: i = 1 estrain@0: n_reads = 0 estrain@0: total_coverage = 0 estrain@0: target_mers = {} estrain@0: for line in open(file): estrain@0: if i % 4 == 2: estrain@0: for kmer, rc_kmer in kmers(line.strip().upper(), k): estrain@0: if (kmer in kmerDict) or (rc_kmer in kmerDict): estrain@0: if kmer in target_mers: estrain@0: target_mers[kmer] += 1 estrain@0: else: estrain@0: target_mers[kmer] = 1 estrain@0: if rc_kmer in target_mers: estrain@0: target_mers[rc_kmer] += 1 estrain@0: else: estrain@0: target_mers[rc_kmer] = 1 estrain@0: i += 1 estrain@0: return set([h for h in target_mers]) estrain@0: estrain@0: estrain@0: def multifasta_single_string2(multifasta): estrain@0: single_string = '' estrain@0: with open(multifasta, 'r') as f: estrain@0: for line in f: estrain@0: if line.strip()[0] == '>': estrain@0: pass estrain@0: else: estrain@0: single_string += line.strip() estrain@0: return single_string estrain@0: estrain@0: estrain@0: def kmers(seq, k): estrain@0: rev_comp = reverse_complement(seq) estrain@0: for start in range(1, len(seq) - k + 1): estrain@0: yield seq[start:start + k], rev_comp[-(start + k):-start] estrain@0: estrain@0: estrain@0: def multifasta_to_kmers_dict(multifasta,k_size):#used to create database kmer set estrain@0: multi_seq_dict = multifasta_dict(multifasta) estrain@0: lib_dict = {} estrain@0: for h in multi_seq_dict: estrain@0: lib_dict[h] = set( estrain@0: [k for k in createKmerDict_reads([multi_seq_dict[h]], k_size)]) estrain@0: return lib_dict estrain@0: estrain@0: estrain@0: def Combine(b, c): estrain@0: fliC_combinations = [] estrain@0: fliC_combinations.append(",".join(c)) estrain@0: temp_combinations = [] estrain@0: for i in range(len(b)): estrain@0: for x in itertools.combinations(b, i + 1): estrain@0: temp_combinations.append(",".join(x)) estrain@0: for x in temp_combinations: estrain@0: temp = [] estrain@0: for y in c: estrain@0: temp.append(y) estrain@0: temp.append(x) estrain@0: temp = ",".join(temp) estrain@0: temp = temp.split(",") estrain@0: temp.sort() estrain@0: temp = ",".join(temp) estrain@0: fliC_combinations.append(temp) estrain@0: return fliC_combinations estrain@0: estrain@0: estrain@0: def seqsero_from_formula_to_serotypes(Otype, fliC, fljB, special_gene_list,subspecies): estrain@0: #like test_output_06012017.txt estrain@0: #can add more varialbles like sdf-type, sub-species-type in future (we can conclude it into a special-gene-list) estrain@0: from Initial_Conditions import phase1,phase2,phaseO,sero,subs,remove_list,rename_dict estrain@0: rename_dict_not_anymore=[rename_dict[x] for x in rename_dict] estrain@0: rename_dict_all=rename_dict_not_anymore+list(rename_dict) #used for decide whether to estrain@0: seronames = [] estrain@0: seronames_none_subspecies=[] estrain@0: for i in range(len(phase1)): estrain@0: fliC_combine = [] estrain@0: fljB_combine = [] estrain@0: if phaseO[i] == Otype: # no VII in KW, but it's there estrain@0: ### for fliC, detect every possible combinations to avoid the effect of "[" estrain@0: if phase1[i].count("[") == 0: estrain@0: fliC_combine.append(phase1[i]) estrain@0: elif phase1[i].count("[") >= 1: estrain@0: c = [] estrain@0: b = [] estrain@0: if phase1[i][0] == "[" and phase1[i][-1] == "]" and phase1[i].count( estrain@0: "[") == 1: estrain@0: content = phase1[i].replace("[", "").replace("]", "") estrain@0: fliC_combine.append(content) estrain@0: fliC_combine.append("-") estrain@0: else: estrain@0: for x in phase1[i].split(","): estrain@0: if "[" in x: estrain@0: b.append(x.replace("[", "").replace("]", "")) estrain@0: else: estrain@0: c.append(x) estrain@0: fliC_combine = Combine( estrain@0: b, c estrain@0: ) #Combine will offer every possible combinations of the formula, like f,[g],t: f,t f,g,t estrain@0: ### end of fliC "[" detect estrain@0: ### for fljB, detect every possible combinations to avoid the effect of "[" estrain@0: if phase2[i].count("[") == 0: estrain@0: fljB_combine.append(phase2[i]) estrain@0: elif phase2[i].count("[") >= 1: estrain@0: d = [] estrain@0: e = [] estrain@0: if phase2[i][0] == "[" and phase2[i][-1] == "]" and phase2[i].count( estrain@0: "[") == 1: estrain@0: content = phase2[i].replace("[", "").replace("]", "") estrain@0: fljB_combine.append(content) estrain@0: fljB_combine.append("-") estrain@0: else: estrain@0: for x in phase2[i].split(","): estrain@0: if "[" in x: estrain@0: d.append(x.replace("[", "").replace("]", "")) estrain@0: else: estrain@0: e.append(x) estrain@0: fljB_combine = Combine(d, e) estrain@0: ### end of fljB "[" detect estrain@0: new_fliC = fliC.split( estrain@0: "," estrain@0: ) #because some antigen like r,[i] not follow alphabetical order, so use this one to judge and can avoid missings estrain@0: new_fliC.sort() estrain@0: new_fliC = ",".join(new_fliC) estrain@0: new_fljB = fljB.split(",") estrain@0: new_fljB.sort() estrain@0: new_fljB = ",".join(new_fljB) estrain@0: if (new_fliC in fliC_combine estrain@0: or fliC in fliC_combine) and (new_fljB in fljB_combine estrain@0: or fljB in fljB_combine): estrain@0: ######start, remove_list,rename_dict, added on 11/11/2018 estrain@0: if sero[i] not in remove_list: estrain@0: temp_sero=sero[i] estrain@0: if temp_sero in rename_dict: estrain@0: temp_sero=rename_dict[temp_sero] #rename if in the rename list estrain@0: if temp_sero not in seronames:#the new sero may already included, if yes, then not consider estrain@0: if subs[i] == subspecies: estrain@0: seronames.append(temp_sero) estrain@0: seronames_none_subspecies.append(temp_sero) estrain@0: else: estrain@0: pass estrain@0: else: estrain@0: pass estrain@0: ######end, added on 11/11/2018 estrain@0: #analyze seronames estrain@0: subspecies_pointer="" estrain@0: if len(seronames) == 0 and len(seronames_none_subspecies)!=0: estrain@0: seronames=seronames_none_subspecies estrain@0: subspecies_pointer="1" estrain@0: if len(seronames) == 0: estrain@0: seronames = [ estrain@0: "N/A (The predicted antigenic profile does not exist in the White-Kauffmann-Le Minor scheme)" estrain@0: ] estrain@0: star = "" estrain@0: star_line = "" estrain@0: if len(seronames) > 1: #there are two possible predictions for serotypes estrain@0: star = "*" estrain@0: #changed 04072019 estrain@0: #star_line = "The predicted serotypes share the same general formula:\t" + Otype + ":" + fliC + ":" + fljB + "\n" estrain@0: if subspecies_pointer=="1" and len(seronames_none_subspecies)!=0: estrain@0: star="*" estrain@0: star_line=" The predicted O and H antigens correspond to serotype '"+(" or ").join(seronames)+"' in the Kauffmann-White scheme. The predicted subspecies by SalmID (github.com/hcdenbakker/SalmID) may not be consistent with subspecies designation in the Kauffmann-White scheme." + star_line estrain@0: #star_line="The formula with this subspieces prediction can't get a serotype in KW manual, and the serotyping prediction was made without considering it."+star_line estrain@0: if Otype=="": estrain@0: Otype="-" estrain@0: predict_form = Otype + ":" + fliC + ":" + fljB estrain@0: predict_sero = (" or ").join(seronames) estrain@0: ###special test for Enteritidis estrain@0: if predict_form == "9:g,m:-": estrain@0: sdf = "-" estrain@0: for x in special_gene_list: estrain@0: if x.startswith("sdf"): estrain@0: sdf = "+" estrain@0: #star_line="Detected sdf gene, a marker to differentiate Gallinarum and Enteritidis" estrain@0: star_line=" sdf gene detected." # ed_SL_04152019: new output format estrain@0: #predict_form = predict_form + " Sdf prediction:" + sdf estrain@0: predict_form = predict_form #changed 04072019 estrain@0: if sdf == "-": estrain@0: star = "*" estrain@0: #star_line="Didn't detected sdf gene, a marker to differentiate Gallinarum and Enteritidis" estrain@0: star_line=" sdf gene not detected." # ed_SL_04152019: new output format estrain@0: #changed in 04072019, for new output estrain@0: #star_line = "Additional characterization is necessary to assign a serotype to this strain. Commonly circulating strains of serotype Enteritidis are sdf+, although sdf- strains of serotype Enteritidis are known to exist. Serotype Gallinarum is typically sdf- but should be quite rare. Sdf- strains of serotype Enteritidis and serotype Gallinarum can be differentiated by phenotypic profile or genetic criteria.\n" estrain@0: #predict_sero = "Gallinarum/Enteritidis" #04132019, for new output requirement estrain@0: predict_sero = "Gallinarum or Enteritidis" # ed_SL_04152019: new output format estrain@0: ###end of special test for Enteritidis estrain@0: elif predict_form == "4:i:-": estrain@0: predict_sero = "I 4,[5],12:i:-" # ed_SL_09242019: change serotype name estrain@0: elif predict_form == "4:r:-": estrain@0: predict_sero = "4:r:-" estrain@0: elif predict_form == "4:b:-": # ed_SL_09272019: change for new output format estrain@0: predict_sero = "N/A (4:b:-)" estrain@0: #elif predict_form == "8:e,h:1,2": #removed after official merge of newport and bardo estrain@0: #predict_sero = "Newport" estrain@0: #star = "*" estrain@0: #star_line = "Serotype Bardo shares the same antigenic profile with Newport, but Bardo is exceedingly rare." estrain@0: claim = " The serotype(s) is/are the only serotype(s) with the indicated antigenic profile currently recognized in the Kauffmann White Scheme. New serotypes can emerge and the possibility exists that this antigenic profile may emerge in a different subspecies. Identification of strains to the subspecies level should accompany serotype determination; the same antigenic profile in different subspecies is considered different serotypes.\n" estrain@0: if "N/A" in predict_sero: estrain@0: claim = "" estrain@0: #special test for Typhimurium estrain@0: if "Typhimurium" in predict_sero or predict_form == "4:i:-": estrain@0: normal = 0 estrain@0: mutation = 0 estrain@0: for x in special_gene_list: estrain@0: if "oafA-O-4_full" in x: estrain@0: normal = float(special_gene_list[x]) estrain@0: elif "oafA-O-4_5-" in x: estrain@0: mutation = float(special_gene_list[x]) estrain@0: if normal > mutation: estrain@0: pass estrain@0: elif normal < mutation: estrain@0: #predict_sero = predict_sero.strip() + "(O5-)" estrain@0: predict_sero = predict_sero.strip() #diable special sero for new output requirement, 04132019 estrain@0: star = "*" estrain@0: #star_line = "Detected the deletion of O5-." estrain@0: star_line = " Detected a deletion that causes O5- variant of Typhimurium." # ed_SL_04152019: new output format estrain@0: else: estrain@0: pass estrain@0: #special test for Paratyphi B estrain@0: if "Paratyphi B" in predict_sero or predict_form == "4:b:-": estrain@0: normal = 0 estrain@0: mutation = 0 estrain@0: for x in special_gene_list: estrain@0: if "gntR-family-regulatory-protein_dt-positive" in x: estrain@0: normal = float(special_gene_list[x]) estrain@0: elif "gntR-family-regulatory-protein_dt-negative" in x: estrain@0: mutation = float(special_gene_list[x]) estrain@0: #print(normal,mutation) estrain@0: if normal > mutation: estrain@0: #predict_sero = predict_sero.strip() + "(dt+)" #diable special sero for new output requirement, 04132019 estrain@0: predict_sero = predict_sero.strip()+' var. L(+) tartrate+' if "Paratyphi B" in predict_sero else predict_sero.strip() # ed_SL_04152019: new output format estrain@0: star = "*" estrain@0: #star_line = "Didn't detect the SNP for dt- which means this isolate is a Paratyphi B variant L(+) tartrate(+)." estrain@0: star_line = " The SNP that causes d-Tartrate nonfermentating phenotype of Paratyphi B was not detected. " # ed_SL_04152019: new output format estrain@0: elif normal < mutation: estrain@0: #predict_sero = predict_sero.strip() + "(dt-)" #diable special sero for new output requirement, 04132019 estrain@0: predict_sero = predict_sero.strip() estrain@0: star = "*" estrain@0: #star_line = "Detected the SNP for dt- which means this isolate is a systemic pathovar of Paratyphi B." estrain@0: star_line = " Detected the SNP for d-Tartrate nonfermenting phenotype of Paratyphi B." # ed_SL_04152019: new output format estrain@0: else: estrain@0: star = "*" estrain@0: #star_line = " Failed to detect the SNP for dt-, can't decide it's a Paratyphi B variant L(+) tartrate(+) or not." estrain@0: star_line = " " ## ed_SL_05152019: do not report this situation. estrain@0: #special test for O13,22 and O13,23 estrain@0: if Otype=="13": estrain@0: #ex_dir = os.path.dirname(os.path.realpath(__file__)) estrain@0: ex_dir = os.path.abspath(os.path.join(os.path.dirname(os.path.dirname(__file__)),'seqsero2_db')) # ed_SL_09152019 estrain@0: f = open(ex_dir + '/special.pickle', 'rb') estrain@0: special = pickle.load(f) estrain@0: O22_O23=special['O22_O23'] estrain@0: if predict_sero.split(" or ")[0] in O22_O23[-1] and predict_sero.split(" or ")[0] not in rename_dict_all:#if in rename_dict_all, then it means already merged, no need to analyze estrain@0: O22_score=0 estrain@0: O23_score=0 estrain@0: for x in special_gene_list: estrain@0: if "O:22" in x: estrain@0: O22_score = O22_score+float(special_gene_list[x]) estrain@0: elif "O:23" in x: estrain@0: O23_score = O23_score+float(special_gene_list[x]) estrain@0: #print(O22_score,O23_score) estrain@0: for z in O22_O23[0]: estrain@0: if predict_sero.split(" or ")[0] in z: estrain@0: if O22_score > O23_score: estrain@0: star = "*" estrain@0: #star_line = "Detected O22 specific genes to further differenciate '"+predict_sero+"'." #diabled for new output requirement, 04132019 estrain@0: predict_sero = z[0] estrain@0: elif O22_score < O23_score: estrain@0: star = "*" estrain@0: #star_line = "Detected O23 specific genes to further differenciate '"+predict_sero+"'." #diabled for new output requirement, 04132019 estrain@0: predict_sero = z[1] estrain@0: else: estrain@0: star = "*" estrain@0: #star_line = "Fail to detect O22 and O23 differences." #diabled for new output requirement, 04132019 estrain@0: if " or " in predict_sero: estrain@0: star_line = star_line + " The predicted serotypes share the same general formula:\t" + Otype + ":" + fliC + ":" + fljB + "\n" estrain@0: #special test for O6,8 estrain@0: #merge_O68_list=["Blockley","Bovismorbificans","Hadar","Litchfield","Manhattan","Muenchen"] #remove 11/11/2018, because already in merge list estrain@0: #for x in merge_O68_list: estrain@0: # if x in predict_sero: estrain@0: # predict_sero=x estrain@0: # star="" estrain@0: # star_line="" estrain@0: #special test for Montevideo; most of them are monophasic estrain@0: #if "Montevideo" in predict_sero and "1,2,7" in predict_form: #remove 11/11/2018, because already in merge list estrain@0: #star="*" estrain@0: #star_line="Montevideo is almost always monophasic, having an antigen called for the fljB position may be a result of Salmonella-Salmonella contamination." estrain@0: return predict_form, predict_sero, star, star_line, claim estrain@0: ### End of SeqSero Kmer part estrain@0: estrain@0: ### Begin of SeqSero2 allele prediction and output estrain@0: def xml_parse_score_comparision_seqsero(xmlfile): estrain@0: #used to do seqsero xml analysis estrain@0: from Bio.Blast import NCBIXML estrain@0: handle=open(xmlfile) estrain@0: handle=NCBIXML.parse(handle) estrain@0: handle=list(handle) estrain@0: List=[] estrain@0: List_score=[] estrain@0: List_ids=[] estrain@0: List_query_region=[] estrain@0: for i in range(len(handle)): estrain@0: if len(handle[i].alignments)>0: estrain@0: for j in range(len(handle[i].alignments)): estrain@0: score=0 estrain@0: ids=0 estrain@0: cover_region=set() #fixed problem that repeated calculation leading percentage > 1 estrain@0: List.append(handle[i].query.strip()+"___"+handle[i].alignments[j].hit_def) estrain@0: for z in range(len(handle[i].alignments[j].hsps)): estrain@0: hsp=handle[i].alignments[j].hsps[z] estrain@0: temp=set(range(hsp.query_start,hsp.query_end)) estrain@0: est_bit=(handle[i].ka_params[0]*hsp.score-math.log(handle[i].ka_params[1]))/(math.log(2)) estrain@0: if len(cover_region)==0: estrain@0: cover_region=cover_region|temp estrain@0: fraction=1 estrain@0: else: estrain@0: fraction=1-len(cover_region&temp)/float(len(temp)) estrain@0: cover_region=cover_region|temp estrain@0: if "last" in handle[i].query or "first" in handle[i].query: estrain@0: #score+=hsp.bits*fraction estrain@0: score+=est_bit*fraction estrain@0: ids+=float(hsp.identities)/handle[i].query_length*fraction estrain@0: else: estrain@0: #score+=hsp.bits*fraction estrain@0: score+=est_bit*fraction estrain@0: ids+=float(hsp.identities)/handle[i].query_length*fraction estrain@0: List_score.append(score) estrain@0: List_ids.append(ids) estrain@0: List_query_region.append(cover_region) estrain@0: temp=zip(List,List_score,List_ids,List_query_region) estrain@0: Final_list=sorted(temp, key=lambda d:d[1], reverse = True) estrain@0: return Final_list estrain@0: estrain@0: estrain@0: def Uniq(L,sort_on_fre="none"): #return the uniq list and the count number estrain@0: Old=L estrain@0: L.sort() estrain@0: L = [L[i] for i in range(len(L)) if L[i] not in L[:i]] estrain@0: count=[] estrain@0: for j in range(len(L)): estrain@0: y=0 estrain@0: for x in Old: estrain@0: if L[j]==x: estrain@0: y+=1 estrain@0: count.append(y) estrain@0: if sort_on_fre!="none": estrain@0: d=zip(*sorted(zip(count, L))) estrain@0: L=d[1] estrain@0: count=d[0] estrain@0: return (L,count) estrain@0: estrain@0: def judge_fliC_or_fljB_from_head_tail_for_one_contig(nodes_vs_score_list): estrain@0: #used to predict it's fliC or fljB for one contig, based on tail and head score, but output the score difference,if it is very small, then not reliable, use blast score for whole contig to test estrain@0: #this is mainly used for estrain@0: a=nodes_vs_score_list estrain@0: fliC_score=0 estrain@0: fljB_score=0 estrain@0: for z in a: estrain@0: if "fliC" in z[0]: estrain@0: fliC_score+=z[1] estrain@0: elif "fljB" in z[0]: estrain@0: fljB_score+=z[1] estrain@0: if fliC_score>=fljB_score: estrain@0: role="fliC" estrain@0: else: estrain@0: role="fljB" estrain@0: return (role,abs(fliC_score-fljB_score)) estrain@0: estrain@0: def judge_fliC_or_fljB_from_whole_contig_blast_score_ranking(node_name,Final_list,Final_list_passed): estrain@0: #used to predict contig is fliC or fljB, if the differnce score value on above head_and_tail is less than 10 (quite small) estrain@0: #also used when no head or tail got blasted score for the contig estrain@0: role="" estrain@0: for z in Final_list_passed: estrain@0: if node_name in z[0]: estrain@0: role=z[0].split("_")[0] estrain@0: break estrain@0: return role estrain@0: estrain@0: def fliC_or_fljB_judge_from_head_tail_sequence(nodes_list,tail_head_list,Final_list,Final_list_passed): estrain@0: #nodes_list is the c created by c,d=Uniq(nodes) in below function estrain@0: first_target="" estrain@0: role_list=[] estrain@0: for x in nodes_list: estrain@0: a=[] estrain@0: role="" estrain@0: for y in tail_head_list: estrain@0: if x in y[0]: estrain@0: a.append(y) estrain@0: if len(a)==4: estrain@0: role,diff=judge_fliC_or_fljB_from_head_tail_for_one_contig(a) estrain@0: if diff<20: estrain@0: role=judge_fliC_or_fljB_from_whole_contig_blast_score_ranking(x,Final_list,Final_list_passed) estrain@0: elif len(a)==3: estrain@0: ###however, if the one with highest score is the fewer one, compare their accumulation score estrain@0: role,diff=judge_fliC_or_fljB_from_head_tail_for_one_contig(a) estrain@0: if diff<20: estrain@0: role=judge_fliC_or_fljB_from_whole_contig_blast_score_ranking(x,Final_list,Final_list_passed) estrain@0: ###end of above score comparison estrain@0: elif len(a)==2: estrain@0: #must on same node, if not, then decide with unit blast score, blast-score/length_of_special_sequence(30 or 37) estrain@0: temp=[] estrain@0: for z in a: estrain@0: temp.append(z[0].split("_")[0]) estrain@0: m,n=Uniq(temp)#should only have one choice, but weird situation might occur too estrain@0: if len(m)==1: estrain@0: pass estrain@0: else: estrain@0: pass estrain@0: role,diff=judge_fliC_or_fljB_from_head_tail_for_one_contig(a) estrain@0: if diff<20: estrain@0: role=judge_fliC_or_fljB_from_whole_contig_blast_score_ranking(x,Final_list,Final_list_passed) estrain@0: ###need to desgin a algorithm to guess most possible situation for nodes_list, See the situations of test evaluation estrain@0: elif len(a)==1: estrain@0: #that one estrain@0: role,diff=judge_fliC_or_fljB_from_head_tail_for_one_contig(a) estrain@0: if diff<20: estrain@0: role=judge_fliC_or_fljB_from_whole_contig_blast_score_ranking(x,Final_list,Final_list_passed) estrain@0: #need to evaluate, in future, may set up a cut-off, if not met, then just find Final_list_passed best match,like when "a==0" estrain@0: else:#a==0 estrain@0: #use Final_list_passed best match estrain@0: for z in Final_list_passed: estrain@0: if x in z[0]: estrain@0: role=z[0].split("_")[0] estrain@0: break estrain@0: #print x,role,len(a) estrain@0: role_list.append((role,x)) estrain@0: if len(role_list)==2: estrain@0: if role_list[0][0]==role_list[1][0]:#this is the most cocmmon error, two antigen were assigned to same phase estrain@0: #just use score to do a final test estrain@0: role_list=[] estrain@0: for x in nodes_list: estrain@0: role=judge_fliC_or_fljB_from_whole_contig_blast_score_ranking(x,Final_list,Final_list_passed) estrain@0: role_list.append((role,x)) estrain@0: return role_list estrain@0: estrain@0: def decide_contig_roles_for_H_antigen(Final_list,Final_list_passed): estrain@0: #used to decide which contig is FliC and which one is fljB estrain@0: contigs=[] estrain@0: nodes=[] estrain@0: for x in Final_list_passed: estrain@0: if x[0].startswith("fl") and "last" not in x[0] and "first" not in x[0]: estrain@0: nodes.append(x[0].split("___")[1].strip()) estrain@0: c,d=Uniq(nodes)#c is node_list estrain@0: #print c estrain@0: tail_head_list=[x for x in Final_list if ("last" in x[0] or "first" in x[0])] estrain@0: roles=fliC_or_fljB_judge_from_head_tail_sequence(c,tail_head_list,Final_list,Final_list_passed) estrain@0: return roles estrain@0: estrain@0: def decide_O_type_and_get_special_genes(Final_list,Final_list_passed): estrain@0: #decide O based on Final_list estrain@0: O_choice="?" estrain@0: O_list=[] estrain@0: special_genes={} estrain@0: nodes=[] estrain@0: for x in Final_list_passed: estrain@0: if x[0].startswith("O-"): estrain@0: nodes.append(x[0].split("___")[1].strip()) estrain@0: elif not x[0].startswith("fl"): estrain@0: special_genes[x[0]]=x[2]#08172018, x[2] changed from x[-1] estrain@0: #print "special_genes:",special_genes estrain@0: c,d=Uniq(nodes) estrain@0: #print "potential O antigen contig",c estrain@0: final_O=[] estrain@0: O_nodes_list=[] estrain@0: for x in c:#c is the list for contigs estrain@0: temp=0 estrain@0: for y in Final_list_passed: estrain@0: if x in y[0] and y[0].startswith("O-"): estrain@0: final_O.append(y) estrain@0: break estrain@0: ### O contig has the problem of two genes on same contig, so do additional test estrain@0: potenial_new_gene="" estrain@0: for x in final_O: estrain@0: pointer=0 #for genes merged or not estrain@0: #not consider O-1,3,19_not_in_3,10, too short compared with others estrain@0: if "O-1,3,19_not_in_3,10" not in x[0] and int(x[0].split("__")[1].split("___")[0])*x[2]+850 <= int(x[0].split("length_")[1].split("_")[0]):#gene length << contig length; for now give 300*2 (for secureity can use 400*2) as flank region estrain@0: pointer=x[0].split("___")[1].strip()#store the contig name estrain@0: print(pointer) estrain@0: if pointer!=0:#it has potential merge event estrain@0: for y in Final_list: estrain@0: if pointer in y[0] and y not in final_O and (y[1]>=int(y[0].split("__")[1].split("___")[0])*1.5 or (y[1]>=int(y[0].split("__")[1].split("___")[0])*y[2] and y[1]>=400)):#that's a realtively strict filter now; if passed, it has merge event and add one more to final_O estrain@0: potenial_new_gene=y estrain@0: #print(potenial_new_gene) estrain@0: break estrain@0: if potenial_new_gene!="": estrain@0: print("two differnt genes in same contig, fix it for O antigen") estrain@0: print(potenial_new_gene[:3]) estrain@0: pointer=0 estrain@0: for y in final_O: estrain@0: if y[0].split("___")[-1]==potenial_new_gene[0].split("___")[-1]: estrain@0: pointer=1 estrain@0: if pointer!=0: #changed to consider two genes in same contig estrain@0: final_O.append(potenial_new_gene) estrain@0: ### end of the two genes on same contig test estrain@0: final_O=sorted(final_O,key=lambda x: x[2], reverse=True)#sorted estrain@0: if len(final_O)==0 or (len(final_O)==1 and "O-1,3,19_not_in_3,10" in final_O[0][0]): estrain@0: #print "$$$No Otype, due to no hit"#may need to be changed estrain@0: O_choice="-" estrain@0: else: estrain@0: highest_O_coverage=max([float(x[0].split("_cov_")[-1].split("_")[0]) for x in final_O if "O-1,3,19_not_in_3,10" not in x[0]]) estrain@0: O_list=[] estrain@0: O_list_less_contamination=[] estrain@0: for x in final_O: estrain@0: if not "O-1,3,19_not_in_3,10__130" in x[0]:#O-1,3,19_not_in_3,10 is too small, which may affect further analysis; to avoid contamination affect, use 0.15 of highest coverage as cut-off estrain@0: O_list.append(x[0].split("__")[0]) estrain@0: O_nodes_list.append(x[0].split("___")[1]) estrain@0: if float(x[0].split("_cov_")[-1].split("_")[0])>highest_O_coverage*0.15: estrain@0: O_list_less_contamination.append(x[0].split("__")[0]) estrain@0: ### special test for O9,46 and O3,10 family estrain@0: if ("O-9,46_wbaV" in O_list or "O-9,46_wbaV-from-II-9,12:z29:1,5-SRR1346254" in O_list) and O_list_less_contamination[0].startswith("O-9,"):#not sure should use and float(O9_wbaV)/float(num_1) > 0.1 estrain@0: if "O-9,46_wzy" in O_list:#and float(O946_wzy)/float(num_1) > 0.1 estrain@0: O_choice="O-9,46" estrain@0: #print "$$$Most possilble Otype: O-9,46" estrain@0: elif "O-9,46,27_partial_wzy" in O_list:#and float(O94627)/float(num_1) > 0.1 estrain@0: O_choice="O-9,46,27" estrain@0: #print "$$$Most possilble Otype: O-9,46,27" estrain@0: else: estrain@0: O_choice="O-9"#next, detect O9 vs O2? estrain@0: O2=0 estrain@0: O9=0 estrain@0: for z in special_genes: estrain@0: if "tyr-O-9" in z: estrain@0: O9=special_genes[z] estrain@0: elif "tyr-O-2" in z: estrain@0: O2=special_genes[z] estrain@0: if O2>O9: estrain@0: O_choice="O-2" estrain@0: elif O2 0.1 and float(O946_wzy)/float(num_1) > 0.1 estrain@0: if "O-3,10_not_in_1,3,19" in O_list:#and float(O310_no_1319)/float(num_1) > 0.1 estrain@0: O_choice="O-3,10" estrain@0: #print "$$$Most possilble Otype: O-3,10 (contain O-3,10_not_in_1,3,19)" estrain@0: else: estrain@0: O_choice="O-1,3,19" estrain@0: #print "$$$Most possilble Otype: O-1,3,19 (not contain O-3,10_not_in_1,3,19)" estrain@0: ### end of special test for O9,46 and O3,10 family estrain@0: else: estrain@0: try: estrain@0: max_score=0 estrain@0: for x in final_O: estrain@0: if x[2]>=max_score and float(x[0].split("_cov_")[-1].split("_")[0])>highest_O_coverage*0.15:#use x[2],08172018, the "coverage identity = cover_length * identity"; also meet coverage threshold estrain@0: max_score=x[2]#change from x[-1] to x[2],08172018 estrain@0: O_choice=x[0].split("_")[0] estrain@0: if O_choice=="O-1,3,19": estrain@0: O_choice=final_O[1][0].split("_")[0] estrain@0: #print "$$$Most possilble Otype: ",O_choice estrain@0: except: estrain@0: pass estrain@0: #print "$$$No suitable Otype, or failure of mapping (please check the quality of raw reads)" estrain@0: #print "O:",O_choice,O_nodes_list estrain@0: Otypes=[] estrain@0: for x in O_list: estrain@0: if x!="O-1,3,19_not_in_3,10": estrain@0: if "O-9,46_" not in x: estrain@0: Otypes.append(x.split("_")[0]) estrain@0: else: estrain@0: Otypes.append(x.split("-from")[0])#O-9,46_wbaV-from-II-9,12:z29:1,5-SRR1346254 estrain@0: #Otypes=[x.split("_")[0] for x in O_list if x!="O-1,3,19_not_in_3,10"] estrain@0: Otypes_uniq,Otypes_fre=Uniq(Otypes) estrain@0: contamination_O="" estrain@0: if O_choice=="O-9,46,27" or O_choice=="O-3,10" or O_choice=="O-1,3,19": estrain@0: if len(Otypes_uniq)>2: estrain@0: contamination_O="potential contamination from O antigen signals" estrain@0: else: estrain@0: if len(Otypes_uniq)>1: estrain@0: if O_choice=="O-4" and len(Otypes_uniq)==2 and "O-9,46,27" in Otypes_uniq: #for special 4,12,27 case such as Bredeney and Schwarzengrund estrain@0: contamination_O="" estrain@0: elif O_choice=="O-9,46" and len(Otypes_uniq)==2 and "O-9,46_wbaV" in Otypes_uniq and "O-9,46_wzy" in Otypes_uniq: #for special 4,12,27 case such as Bredeney and Schwarzengrund estrain@0: contamination_O="" estrain@0: else: estrain@0: contamination_O="potential contamination from O antigen signals" estrain@0: return O_choice,O_nodes_list,special_genes,final_O,contamination_O,Otypes_uniq estrain@0: ### End of SeqSero2 allele prediction and output estrain@0: estrain@0: def get_input_files(make_dir,input_file,data_type,dirpath): estrain@0: #tell input files from datatype estrain@0: #": '1'(pair-end reads, interleaved),'2'(pair-end reads, seperated),'3'(single-end reads), '4'(assembly),'5'(nanopore fasta),'6'(nanopore fastq)" estrain@0: for_fq="" estrain@0: rev_fq="" estrain@0: os.chdir(make_dir) estrain@0: if data_type=="1": estrain@0: input_file=input_file[0].split("/")[-1] estrain@0: if input_file.endswith(".sra"): estrain@0: subprocess.check_call("fastq-dump --split-files "+input_file,shell=True) estrain@0: for_fq=input_file.replace(".sra","_1.fastq") estrain@0: rev_fq=input_file.replace(".sra","_2.fastq") estrain@0: else: estrain@0: core_id=input_file.split(".fastq")[0].split(".fq")[0] estrain@0: for_fq=core_id+"_1.fastq" estrain@0: rev_fq=core_id+"_2.fastq" estrain@0: if input_file.endswith(".gz"): estrain@0: subprocess.check_call("gzip -dc "+input_file+" | "+dirpath+"/deinterleave_fastq.sh "+for_fq+" "+rev_fq,shell=True) estrain@0: else: estrain@0: subprocess.check_call("cat "+input_file+" | "+dirpath+"/deinterleave_fastq.sh "+for_fq+" "+rev_fq,shell=True) estrain@0: elif data_type=="2": estrain@0: for_fq=input_file[0].split("/")[-1] estrain@0: rev_fq=input_file[1].split("/")[-1] estrain@0: elif data_type=="3": estrain@0: input_file=input_file[0].split("/")[-1] estrain@0: if input_file.endswith(".sra"): estrain@0: subprocess.check_call("fastq-dump --split-files "+input_file,shell=True) estrain@0: for_fq=input_file.replace(".sra","_1.fastq") estrain@0: else: estrain@0: for_fq=input_file estrain@0: elif data_type in ["4","5","6"]: estrain@0: for_fq=input_file[0].split("/")[-1] estrain@0: os.chdir("..") estrain@0: return for_fq,rev_fq estrain@0: estrain@0: def predict_O_and_H_types(Final_list,Final_list_passed,new_fasta): estrain@0: #get O and H types from Final_list from blast parsing; allele mode estrain@0: from Bio import SeqIO estrain@0: fliC_choice="-" estrain@0: fljB_choice="-" estrain@0: fliC_contig="NA" estrain@0: fljB_contig="NA" estrain@0: fliC_region=set([0]) estrain@0: fljB_region=set([0,]) estrain@0: fliC_length=0 #can be changed to coverage in future; in 03292019, changed to ailgned length estrain@0: fljB_length=0 #can be changed to coverage in future; in 03292019, changed to ailgned length estrain@0: O_choice="-"#no need to decide O contig for now, should be only one estrain@0: O_choice,O_nodes,special_gene_list,O_nodes_roles,contamination_O,Otypes_uniq=decide_O_type_and_get_special_genes(Final_list,Final_list_passed)#decide the O antigen type and also return special-gene-list for further identification estrain@0: O_choice=O_choice.split("-")[-1].strip() estrain@0: if (O_choice=="1,3,19" and len(O_nodes_roles)==1 and "1,3,19" in O_nodes_roles[0][0]) or O_choice=="": estrain@0: O_choice="-" estrain@0: H_contig_roles=decide_contig_roles_for_H_antigen(Final_list,Final_list_passed)#decide the H antigen contig is fliC or fljB estrain@0: #add alignment locations, used for further selection, 03312019 estrain@0: for i in range(len(H_contig_roles)): estrain@0: x=H_contig_roles[i] estrain@0: for y in Final_list_passed: estrain@0: if x[1] in y[0] and y[0].startswith(x[0]): estrain@0: H_contig_roles[i]+=H_contig_roles[i]+(y[-1],) estrain@0: break estrain@0: log_file=open("SeqSero_log.txt","a") estrain@0: extract_file=open("Extracted_antigen_alleles.fasta","a") estrain@0: handle_fasta=list(SeqIO.parse(new_fasta,"fasta")) estrain@0: estrain@0: #print("O_contigs:") estrain@0: log_file.write("O_contigs:\n") estrain@0: extract_file.write("#Sequences with antigen signals (if the micro-assembled contig only covers the flanking region, it will not be used for contamination analysis)\n") estrain@0: extract_file.write("#O_contigs:\n") estrain@0: for x in O_nodes_roles: estrain@0: if "O-1,3,19_not_in_3,10" not in x[0]:#O-1,3,19_not_in_3,10 is just a small size marker estrain@0: #print(x[0].split("___")[-1],x[0].split("__")[0],"blast score:",x[1],"identity%:",str(round(x[2]*100,2))+"%",str(min(x[-1]))+" to "+str(max(x[-1]))) estrain@0: log_file.write(x[0].split("___")[-1]+" "+x[0].split("__")[0]+"; "+"blast score: "+str(x[1])+" identity%: "+str(round(x[2]*100,2))+"%; alignment from "+str(min(x[-1]))+" to "+str(max(x[-1]))+" of antigen\n") estrain@0: title=">"+x[0].split("___")[-1]+" "+x[0].split("__")[0]+"; "+"blast score: "+str(x[1])+" identity%: "+str(round(x[2]*100,2))+"%; alignment from "+str(min(x[-1]))+" to "+str(max(x[-1]))+" of antigen\n" estrain@0: seqs="" estrain@0: for z in handle_fasta: estrain@0: if x[0].split("___")[-1]==z.description: estrain@0: seqs=str(z.seq) estrain@0: extract_file.write(title+seqs+"\n") estrain@0: if len(H_contig_roles)!=0: estrain@0: highest_H_coverage=max([float(x[1].split("_cov_")[-1].split("_")[0]) for x in H_contig_roles]) #less than highest*0.1 would be regarded as contamination and noises, they will still be considered in contamination detection and logs, but not used as final serotype output estrain@0: else: estrain@0: highest_H_coverage=0 estrain@0: for x in H_contig_roles: estrain@0: #if multiple choices, temporately select the one with longest length for now, will revise in further change estrain@0: if "fliC" == x[0] and len(x[-1])>=fliC_length and x[1] not in O_nodes and float(x[1].split("_cov_")[-1].split("_")[0])>highest_H_coverage*0.13:#remember to avoid the effect of O-type contig, so should not in O_node list estrain@0: fliC_contig=x[1] estrain@0: fliC_length=len(x[-1]) estrain@0: elif "fljB" == x[0] and len(x[-1])>=fljB_length and x[1] not in O_nodes and float(x[1].split("_cov_")[-1].split("_")[0])>highest_H_coverage*0.13: estrain@0: fljB_contig=x[1] estrain@0: fljB_length=len(x[-1]) estrain@0: for x in Final_list_passed: estrain@0: if fliC_choice=="-" and "fliC_" in x[0] and fliC_contig in x[0]: estrain@0: fliC_choice=x[0].split("_")[1] estrain@0: elif fljB_choice=="-" and "fljB_" in x[0] and fljB_contig in x[0]: estrain@0: fljB_choice=x[0].split("_")[1] estrain@0: elif fliC_choice!="-" and fljB_choice!="-": estrain@0: break estrain@0: #now remove contigs not in middle core part estrain@0: first_allele="NA" estrain@0: first_allele_percentage=0 estrain@0: for x in Final_list: estrain@0: if x[0].startswith("fliC") or x[0].startswith("fljB"): estrain@0: first_allele=x[0].split("__")[0] #used to filter those un-middle contigs estrain@0: first_allele_percentage=x[2] estrain@0: break estrain@0: additional_contigs=[] estrain@0: for x in Final_list: estrain@0: if first_allele in x[0]: estrain@0: if (fliC_contig == x[0].split("___")[-1]): estrain@0: fliC_region=x[3] estrain@0: elif fljB_contig!="NA" and (fljB_contig == x[0].split("___")[-1]): estrain@0: fljB_region=x[3] estrain@0: else: estrain@0: if x[1]*1.1>int(x[0].split("___")[1].split("_")[3]):#loose threshold by multiplying 1.1 estrain@0: additional_contigs.append(x) estrain@0: #else: estrain@0: #print x[:3] estrain@0: #we can just use the fljB region (or fliC depends on size), no matter set() or contain a large locations (without middle part); however, if none of them is fully assembled, use 500 and 1200 as conservative cut-off estrain@0: if first_allele_percentage>0.9: estrain@0: if len(fliC_region)>len(fljB_region) and (max(fljB_region)-min(fljB_region))>1000: estrain@0: target_region=fljB_region|(fliC_region-set(range(min(fljB_region),max(fljB_region)))) #fljB_region|(fliC_region-set(range(min(fljB_region),max(fljB_region)))) estrain@0: elif len(fliC_region)1000: estrain@0: target_region=fliC_region|(fljB_region-set(range(min(fliC_region),max(fliC_region)))) #fljB_region|(fliC_region-set(range(min(fljB_region),max(fljB_region)))) estrain@0: else: estrain@0: target_region=set()#doesn't do anything estrain@0: else: estrain@0: target_region=set()#doesn't do anything estrain@0: #print(target_region) estrain@0: #print(additional_contigs) estrain@0: target_region2=set(list(range(0,525))+list(range(1200,1700)))#I found to use 500 to 1200 as special region would be best estrain@0: target_region=target_region2|target_region estrain@0: for x in additional_contigs: estrain@0: removal=0 estrain@0: contig_length=int(x[0].split("___")[1].split("length_")[-1].split("_")[0]) estrain@0: if fljB_contig not in x[0] and fliC_contig not in x[0] and len(target_region&x[3])/float(len(x[3]))>0.65 and contig_length*0.5 0.9 and float(x[0].split("__")[1].split("___")[0])*x[2]/len(x[-1])>0.96:#if high similiarity with middle part of first allele (first allele >0.9, already cover middle part) estrain@0: removal=1 estrain@0: else: estrain@0: pass estrain@0: if removal==1: estrain@0: for y in H_contig_roles: estrain@0: if y[1] in x[0]: estrain@0: H_contig_roles.remove(y) estrain@0: else: estrain@0: pass estrain@0: #print(x[:3],contig_length,len(target_region&x[3])/float(len(x[3])),contig_length*0.5,len(x[3]),contig_length*1.5) estrain@0: #end of removing none-middle contigs estrain@0: #print("H_contigs:") estrain@0: log_file.write("H_contigs:\n") estrain@0: extract_file.write("#H_contigs:\n") estrain@0: H_contig_stat=[] estrain@0: H1_cont_stat={} estrain@0: H2_cont_stat={} estrain@0: for i in range(len(H_contig_roles)): estrain@0: x=H_contig_roles[i] estrain@0: a=0 estrain@0: for y in Final_list_passed: estrain@0: if x[1] in y[0] and y[0].startswith(x[0]): estrain@0: if "first" in y[0] or "last" in y[0]: #this is the final filter to decide it's fliC or fljB, if can't pass, then can't decide estrain@0: for y in Final_list_passed: #it's impossible to has the "first" and "last" allele as prediction, so re-do it estrain@0: if x[1] in y[0]:#it's very possible to be third phase allele, so no need to make it must be fliC or fljB estrain@0: #print(x[1],"can't_decide_fliC_or_fljB",y[0].split("_")[1],"blast_score:",y[1],"identity%:",str(round(y[2]*100,2))+"%",str(min(y[-1]))+" to "+str(max(y[-1]))) estrain@0: log_file.write(x[1]+" "+x[0]+" "+y[0].split("_")[1]+"; "+"blast score: "+str(y[1])+" identity%: "+str(round(y[2]*100,2))+"%; alignment from "+str(min(y[-1]))+" to "+str(max(y[-1]))+" of antigen\n") estrain@0: H_contig_roles[i]="can't decide fliC or fljB, may be third phase" estrain@0: title=">"+x[1]+" "+x[0]+" "+y[0].split("_")[1]+"; "+"blast score: "+str(y[1])+" identity%: "+str(round(y[2]*100,2))+"%; alignment from "+str(min(y[-1]))+" to "+str(max(y[-1]))+" of antiten\n" estrain@0: seqs="" estrain@0: for z in handle_fasta: estrain@0: if x[1]==z.description: estrain@0: seqs=str(z.seq) estrain@0: extract_file.write(title+seqs+"\n") estrain@0: break estrain@0: else: estrain@0: #print(x[1],x[0],y[0].split("_")[1],"blast_score:",y[1],"identity%:",str(round(y[2]*100,2))+"%",str(min(y[-1]))+" to "+str(max(y[-1]))) estrain@0: log_file.write(x[1]+" "+x[0]+" "+y[0].split("_")[1]+"; "+"blast score: "+str(y[1])+" identity%: "+str(round(y[2]*100,2))+"%; alignment from "+str(min(y[-1]))+" to "+str(max(y[-1]))+" of antigen\n") estrain@0: title=">"+x[1]+" "+x[0]+" "+y[0].split("_")[1]+"; "+"blast score: "+str(y[1])+" identity%: "+str(round(y[2]*100,2))+"%; alignment from "+str(min(y[-1]))+" to "+str(max(y[-1]))+" of antigen\n" estrain@0: seqs="" estrain@0: for z in handle_fasta: estrain@0: if x[1]==z.description: estrain@0: seqs=str(z.seq) estrain@0: extract_file.write(title+seqs+"\n") estrain@0: if x[0]=="fliC": estrain@0: if y[0].split("_")[1] not in H1_cont_stat: estrain@0: H1_cont_stat[y[0].split("_")[1]]=y[2] estrain@0: else: estrain@0: H1_cont_stat[y[0].split("_")[1]]+=y[2] estrain@0: if x[0]=="fljB": estrain@0: if y[0].split("_")[1] not in H2_cont_stat: estrain@0: H2_cont_stat[y[0].split("_")[1]]=y[2] estrain@0: else: estrain@0: H2_cont_stat[y[0].split("_")[1]]+=y[2] estrain@0: break estrain@0: #detect contaminations estrain@0: #print(H1_cont_stat) estrain@0: #print(H2_cont_stat) estrain@0: H1_cont_stat_list=[x for x in H1_cont_stat if H1_cont_stat[x]>0.2] estrain@0: H2_cont_stat_list=[x for x in H2_cont_stat if H2_cont_stat[x]>0.2] estrain@0: contamination_H="" estrain@0: if len(H1_cont_stat_list)>1 or len(H2_cont_stat_list)>1: estrain@0: contamination_H="potential contamination from H antigen signals" estrain@0: elif len(H2_cont_stat_list)==1 and fljB_contig=="NA": estrain@0: contamination_H="potential contamination from H antigen signals, uncommon weak fljB signals detected" estrain@0: #get additional antigens estrain@0: """ estrain@0: if ("O-9,46_wbaV" in O_list or "O-9,46_wbaV-from-II-9,12:z29:1,5-SRR1346254" in O_list) and O_list_less_contamination[0].startswith("O-9,"):#not sure should use and float(O9_wbaV)/float(num_1) > 0.1 estrain@0: if "O-9,46_wzy" in O_list:#and float(O946_wzy)/float(num_1) > 0.1 estrain@0: O_choice="O-9,46" estrain@0: #print "$$$Most possilble Otype: O-9,46" estrain@0: elif "O-9,46,27_partial_wzy" in O_list:#and float(O94627)/float(num_1) > 0.1 estrain@0: O_choice="O-9,46,27" estrain@0: #print "$$$Most possilble Otype: O-9,46,27" estrain@0: elif ("O-3,10_wzx" in O_list) and ("O-9,46_wzy" in O_list) and (O_list[0].startswith("O-3,10") or O_list_less_contamination[0].startswith("O-9,46_wzy")):#and float(O310_wzx)/float(num_1) > 0.1 and float(O946_wzy)/float(num_1) > 0.1 estrain@0: if "O-3,10_not_in_1,3,19" in O_list:#and float(O310_no_1319)/float(num_1) > 0.1 estrain@0: O_choice="O-3,10" estrain@0: #print "$$$Most possilble Otype: O-3,10 (contain O-3,10_not_in_1,3,19)" estrain@0: else: estrain@0: O_choice="O-1,3,19" estrain@0: #print "$$$Most possilble Otype: O-1,3,19 (not contain O-3,10_not_in_1,3,19)" estrain@0: ### end of special test for O9,46 and O3,10 family estrain@0: estrain@0: if O_choice=="O-9,46,27" or O_choice=="O-3,10" or O_choice=="O-1,3,19": estrain@0: if len(Otypes_uniq)>2: estrain@0: contamination_O="potential contamination from O antigen signals" estrain@0: else: estrain@0: if len(Otypes_uniq)>1: estrain@0: if O_choice=="O-4" and len(Otypes_uniq)==2 and "O-9,46,27" in Otypes_uniq: #for special 4,12,27 case such as Bredeney and Schwarzengrund estrain@0: contamination_O="" estrain@0: elif O_choice=="O-9,46" and len(Otypes_uniq)==2 and "O-9,46_wbaV" in Otypes_uniq and "O-9,46_wzy" in Otypes_uniq: #for special 4,12,27 case such as Bredeney and Schwarzengrund estrain@0: contamination_O="" estrain@0: """ estrain@0: additonal_antigents=[] estrain@0: #print(contamination_O) estrain@0: #print(contamination_H) estrain@0: log_file.write(contamination_O+"\n") estrain@0: log_file.write(contamination_H+"\n") estrain@0: log_file.close() estrain@0: return O_choice,fliC_choice,fljB_choice,special_gene_list,contamination_O,contamination_H,Otypes_uniq,H1_cont_stat_list,H2_cont_stat_list estrain@0: estrain@0: def get_input_K(input_file,lib_dict,data_type,k_size): estrain@0: #kmer mode; get input_Ks from dict and data_type estrain@0: kmers = [] estrain@0: for h in lib_dict: estrain@0: kmers += lib_dict[h] estrain@0: if data_type == '4': estrain@0: input_Ks = target_multifasta_kmerizer(input_file, k_size, set(kmers)) estrain@0: elif data_type == '1' or data_type == '2' or data_type == '3':#set it for now, will change later estrain@0: input_Ks = target_read_kmerizer(input_file, k_size, set(kmers)) estrain@0: elif data_type == '5':#minion_2d_fasta estrain@0: input_Ks = minion_fasta_kmerizer(input_file, k_size, set(kmers)) estrain@0: if data_type == '6':#minion_2d_fastq estrain@0: input_Ks = minion_fastq_kmerizer(input_file, k_size, set(kmers)) estrain@0: return input_Ks estrain@0: estrain@0: def get_kmer_dict(lib_dict,input_Ks): estrain@0: #kmer mode; get predicted types estrain@0: O_dict = {} estrain@0: H_dict = {} estrain@0: Special_dict = {} estrain@0: for h in lib_dict: estrain@0: score = (len(lib_dict[h] & input_Ks) / len(lib_dict[h])) * 100 estrain@0: if score > 1: # Arbitrary cut-off for similarity score very low but seems necessary to detect O-3,10 in some cases estrain@0: if h.startswith('O-') and score > 25: estrain@0: O_dict[h] = score estrain@0: if h.startswith('fl') and score > 40: estrain@0: H_dict[h] = score estrain@0: if (h[:2] != 'fl') and (h[:2] != 'O-'): estrain@0: Special_dict[h] = score estrain@0: return O_dict,H_dict,Special_dict estrain@0: estrain@0: def call_O_and_H_type(O_dict,H_dict,Special_dict,make_dir): estrain@0: log_file=open("SeqSero_log.txt","a") estrain@0: log_file.write("O_scores:\n") estrain@0: #call O: estrain@0: highest_O = '-' estrain@0: if len(O_dict) == 0: estrain@0: pass estrain@0: else: estrain@0: for x in O_dict: estrain@0: log_file.write(x+"\t"+str(O_dict[x])+"\n") estrain@0: if ('O-9,46_wbaV__1002' in O_dict and O_dict['O-9,46_wbaV__1002']>70) or ("O-9,46_wbaV-from-II-9,12:z29:1,5-SRR1346254__1002" in O_dict and O_dict['O-9,46_wbaV-from-II-9,12:z29:1,5-SRR1346254__1002']>70): # not sure should use and float(O9_wbaV)/float(num_1) > 0.1 estrain@0: if 'O-9,46_wzy__1191' in O_dict: # and float(O946_wzy)/float(num_1) > 0.1 estrain@0: highest_O = "O-9,46" estrain@0: elif "O-9,46,27_partial_wzy__1019" in O_dict: # and float(O94627)/float(num_1) > 0.1 estrain@0: highest_O = "O-9,46,27" estrain@0: else: estrain@0: highest_O = "O-9" # next, detect O9 vs O2? estrain@0: O2 = 0 estrain@0: O9 = 0 estrain@0: for z in Special_dict: estrain@0: if "tyr-O-9" in z: estrain@0: O9 = float(Special_dict[z]) estrain@0: if "tyr-O-2" in z: estrain@0: O2 = float(Special_dict[z]) estrain@0: if O2 > O9: estrain@0: highest_O = "O-2" estrain@0: elif ("O-3,10_wzx__1539" in O_dict) and ( estrain@0: "O-9,46_wzy__1191" in O_dict estrain@0: ): # and float(O310_wzx)/float(num_1) > 0.1 and float(O946_wzy)/float(num_1) > 0.1 estrain@0: if "O-3,10_not_in_1,3,19__1519" in O_dict: # and float(O310_no_1319)/float(num_1) > 0.1 estrain@0: highest_O = "O-3,10" estrain@0: else: estrain@0: highest_O = "O-1,3,19" estrain@0: ### end of special test for O9,46 and O3,10 family estrain@0: else: estrain@0: try: estrain@0: max_score = 0 estrain@0: for x in O_dict: estrain@0: if float(O_dict[x]) >= max_score: estrain@0: max_score = float(O_dict[x]) estrain@0: highest_O = x.split("_")[0] estrain@0: if highest_O == "O-1,3,19": estrain@0: highest_O = '-' estrain@0: max_score = 0 estrain@0: for x in O_dict: estrain@0: if x == 'O-1,3,19_not_in_3,10__130': estrain@0: pass estrain@0: else: estrain@0: if float(O_dict[x]) >= max_score: estrain@0: max_score = float(O_dict[x]) estrain@0: highest_O = x.split("_")[0] estrain@0: except: estrain@0: pass estrain@0: #call_fliC: estrain@0: if len(H_dict)!=0: estrain@0: highest_H_score_both_BC=H_dict[max(H_dict.keys(), key=(lambda k: H_dict[k]))] #used to detect whether fljB existed or not estrain@0: else: estrain@0: highest_H_score_both_BC=0 estrain@0: highest_fliC = '-' estrain@0: highest_fliC_raw = '-' estrain@0: highest_Score = 0 estrain@0: log_file.write("\nH_scores:\n") estrain@0: for s in H_dict: estrain@0: log_file.write(s+"\t"+str(H_dict[s])+"\n") estrain@0: if s.startswith('fliC'): estrain@0: if float(H_dict[s]) > highest_Score: estrain@0: highest_fliC = s.split('_')[1] estrain@0: highest_fliC_raw = s estrain@0: highest_Score = float(H_dict[s]) estrain@0: #call_fljB estrain@0: highest_fljB = '-' estrain@0: highest_fljB_raw = '-' estrain@0: highest_Score = 0 estrain@0: for s in H_dict: estrain@0: if s.startswith('fljB'): estrain@0: if float(H_dict[s]) > highest_Score and float(H_dict[s]) > highest_H_score_both_BC * 0.65: #fljB is special, so use highest_H_score_both_BC to give a general estimate of coverage, currently 0.65 seems pretty good; the reason use a high (0.65) is some fliC and fljB shared with each other estrain@0: highest_fljB = s.split('_')[1] estrain@0: highest_fljB_raw = s estrain@0: highest_Score = float(H_dict[s]) estrain@0: log_file.write("\nSpecial_scores:\n") estrain@0: for s in Special_dict: estrain@0: log_file.write(s+"\t"+str(Special_dict[s])+"\n") estrain@0: log_file.close() estrain@0: return highest_O,highest_fliC,highest_fljB estrain@0: estrain@0: def get_temp_file_names(for_fq,rev_fq): estrain@0: #seqsero2 -a; get temp file names estrain@0: sam=for_fq+".sam" estrain@0: bam=for_fq+".bam" estrain@0: sorted_bam=for_fq+"_sorted.bam" estrain@0: mapped_fq1=for_fq+"_mapped.fq" estrain@0: mapped_fq2=rev_fq+"_mapped.fq" estrain@0: combined_fq=for_fq+"_combined.fq" estrain@0: for_sai=for_fq+".sai" estrain@0: rev_sai=rev_fq+".sai" estrain@0: return sam,bam,sorted_bam,mapped_fq1,mapped_fq2,combined_fq,for_sai,rev_sai estrain@0: estrain@0: def map_and_sort(threads,database,fnameA,fnameB,sam,bam,for_sai,rev_sai,sorted_bam,mapping_mode): estrain@0: #seqsero2 -a; do mapping and sort estrain@0: print("building database...") estrain@0: subprocess.check_call("bwa index "+database+ " 2>> data_log.txt",shell=True) estrain@0: print("mapping...") estrain@0: if mapping_mode=="mem": estrain@0: subprocess.check_call("bwa mem -k 17 -t "+threads+" "+database+" "+fnameA+" "+fnameB+" > "+sam+ " 2>> data_log.txt",shell=True) estrain@0: elif mapping_mode=="sam": estrain@0: if fnameB!="": estrain@0: subprocess.check_call("bwa aln -t "+threads+" "+database+" "+fnameA+" > "+for_sai+ " 2>> data_log.txt",shell=True) estrain@0: subprocess.check_call("bwa aln -t "+threads+" "+database+" "+fnameB+" > "+rev_sai+ " 2>> data_log.txt",shell=True) estrain@0: subprocess.check_call("bwa sampe "+database+" "+for_sai+" "+ rev_sai+" "+fnameA+" "+fnameB+" > "+sam+ " 2>> data_log.txt",shell=True) estrain@0: else: estrain@0: subprocess.check_call("bwa aln -t "+threads+" "+database+" "+fnameA+" > "+for_sai+ " 2>> data_log.txt",shell=True) estrain@0: subprocess.check_call("bwa samse "+database+" "+for_sai+" "+for_fq+" > "+sam) estrain@0: subprocess.check_call("samtools view -@ "+threads+" -F 4 -Sh "+sam+" > "+bam,shell=True) estrain@0: ### check the version of samtools then use differnt commands estrain@0: samtools_version=subprocess.Popen(["samtools"],stdout=subprocess.PIPE,stderr=subprocess.PIPE) estrain@0: out, err = samtools_version.communicate() estrain@0: version = str(err).split("ersion:")[1].strip().split(" ")[0].strip() estrain@0: print("check samtools version:",version) estrain@0: ### end of samtools version check and its analysis estrain@0: if LooseVersion(version)<=LooseVersion("1.2"): estrain@0: subprocess.check_call("samtools sort -@ "+threads+" -n "+bam+" "+fnameA+"_sorted",shell=True) estrain@0: else: estrain@0: subprocess.check_call("samtools sort -@ "+threads+" -n "+bam+" >"+sorted_bam,shell=True) estrain@0: estrain@0: def extract_mapped_reads_and_do_assembly_and_blast(current_time,sorted_bam,combined_fq,mapped_fq1,mapped_fq2,threads,fnameA,fnameB,database,mapping_mode): estrain@0: #seqsero2 -a; extract, assembly and blast estrain@0: subprocess.check_call("bamToFastq -i "+sorted_bam+" -fq "+combined_fq,shell=True) estrain@0: #print("fnameA:",fnameA) estrain@0: #print("fnameB:",fnameB) estrain@0: if fnameB!="": estrain@0: subprocess.check_call("bamToFastq -i "+sorted_bam+" -fq "+mapped_fq1+" -fq2 "+mapped_fq2 + " 2>> data_log.txt",shell=True)#2> /dev/null if want no output estrain@0: else: estrain@0: pass estrain@0: outdir=current_time+"_temp" estrain@0: print("assembling...") estrain@0: if int(threads)>4: estrain@0: t="4" estrain@0: else: estrain@0: t=threads estrain@0: if os.path.getsize(combined_fq)>100 and (fnameB=="" or os.path.getsize(mapped_fq1)>100):#if not, then it's "-:-:-" estrain@0: if fnameB!="": estrain@0: subprocess.check_call("spades.py --careful --pe1-s "+combined_fq+" --pe1-1 "+mapped_fq1+" --pe1-2 "+mapped_fq2+" -t "+t+" -o "+outdir+ " >> data_log.txt 2>&1",shell=True) estrain@0: else: estrain@0: subprocess.check_call("spades.py --careful --pe1-s "+combined_fq+" -t "+t+" -o "+outdir+ " >> data_log.txt 2>&1",shell=True) estrain@0: #new_fasta=fnameA+"_"+database+"_"+mapping_mode+".fasta" estrain@0: new_fasta=fnameA+"_"+database.split('/')[-1]+"_"+mapping_mode+".fasta" # ed_SL_09152019: change path to databse for packaging estrain@0: subprocess.check_call("mv "+outdir+"/contigs.fasta "+new_fasta+ " 2> /dev/null",shell=True) estrain@0: #os.system("mv "+outdir+"/scaffolds.fasta "+new_fasta+ " 2> /dev/null") contigs.fasta estrain@0: subprocess.check_call("rm -rf "+outdir+ " 2> /dev/null",shell=True) estrain@0: print("blasting...","\n") estrain@0: xmlfile="blasted_output.xml"#fnameA+"-extracted_vs_"+database+"_"+mapping_mode+".xml" estrain@0: subprocess.check_call('makeblastdb -in '+new_fasta+' -out '+new_fasta+'_db '+'-dbtype nucl >> data_log.txt 2>&1',shell=True) #temp.txt is to forbid the blast result interrupt the output of our program###1/27/2015 estrain@0: subprocess.check_call("blastn -query "+database+" -db "+new_fasta+"_db -out "+xmlfile+" -outfmt 5 >> data_log.txt 2>&1",shell=True)###1/27/2015; 08272018, remove "-word_size 10" estrain@0: else: estrain@0: xmlfile="NA" estrain@0: return xmlfile,new_fasta estrain@0: estrain@0: def judge_subspecies(fnameA,dirpath): estrain@0: #seqsero2 -a; judge subspecies on just forward raw reads fastq estrain@0: salmID_output=subprocess.Popen("python " + dirpath + "/SalmID.py -i "+fnameA,shell=True,stdout=subprocess.PIPE,stderr=subprocess.PIPE) estrain@0: out, err = salmID_output.communicate() estrain@0: out=out.decode("utf-8") estrain@0: file=open("data_log.txt","a") estrain@0: file.write(out) estrain@0: file.close() estrain@0: salm_species_scores=out.split("\n")[1].split("\t")[6:] estrain@0: salm_species_results=out.split("\n")[0].split("\t")[6:] estrain@0: max_score=0 estrain@0: max_score_index=1 #default is 1, means "I" estrain@0: for i in range(len(salm_species_scores)): estrain@0: if max_score float(out.split("\n")[1].split("\t")[5]): #bongori and enterica compare estrain@0: prediction="bongori" #if not, the prediction would always be enterica, since they are located in the later part estrain@0: if max_score<10: estrain@0: prediction="-" estrain@0: return prediction estrain@0: estrain@0: def judge_subspecies_Kmer(Special_dict): estrain@0: #seqsero2 -k; estrain@0: max_score=0 estrain@0: prediction="-" #default should be I estrain@0: for x in Special_dict: estrain@0: if "mer" in x: estrain@0: if max_score95:#if bongori already, then no need to test enterica estrain@0: prediction="bongori" estrain@0: break estrain@0: return prediction estrain@0: estrain@0: def main(): estrain@0: #combine SeqSeroK and SeqSero2, also with SalmID estrain@0: args = parse_args() estrain@0: input_file = args.i estrain@0: data_type = args.t estrain@0: analysis_mode = args.m estrain@0: mapping_mode=args.b estrain@0: threads=args.p estrain@0: make_dir=args.d estrain@0: clean_mode=args.c estrain@0: k_size=27 #will change for bug fixing estrain@0: #database="H_and_O_and_specific_genes.fasta" estrain@0: dirpath = os.path.abspath(os.path.dirname(os.path.realpath(__file__))) estrain@0: ex_dir = os.path.abspath(os.path.join(os.path.dirname(os.path.dirname(__file__)),'seqsero2_db')) # ed_SL_09152019: add ex_dir for packaging estrain@4: database=dirpath+"/seqsero2_db/H_and_O_and_specific_genes.fasta" # ed_SL_09152019: change path to database for packaging estrain@0: note="Note:" estrain@0: NA_note=" This predicted serotype is not in the Kauffman-White scheme." # ed_SL_09272019: add for new output format estrain@0: if len(sys.argv)==1: estrain@0: subprocess.check_call(dirpath+"/SeqSero2_package.py -h",shell=True)#change name of python file estrain@0: else: estrain@0: request_id = time.strftime("%m_%d_%Y_%H_%M_%S", time.localtime()) estrain@0: request_id += str(random.randint(1, 10000000)) estrain@0: if make_dir is None: estrain@0: make_dir="SeqSero_result_"+request_id estrain@0: if os.path.isdir(make_dir): estrain@0: pass estrain@0: else: estrain@0: subprocess.check_call(["mkdir",make_dir]) estrain@0: #subprocess.check_call("cp "+dirpath+"/"+database+" "+" ".join(input_file)+" "+make_dir,shell=True) estrain@0: #subprocess.check_call("ln -sr "+dirpath+"/"+database+" "+" ".join(input_file)+" "+make_dir,shell=True) estrain@0: subprocess.check_call("ln -f -s "+database+" "+" ".join(input_file)+" "+make_dir,shell=True) # ed_SL_09152019: change path to database for packaging estrain@0: #subprocess.check_call("ln -f -s "+dirpath+"/"+database+" "+" ".join(input_file)+" "+make_dir,shell=True) ### ed_SL_05282019: use -f option to force the replacement of links, remove -r and use absolute path instead to avoid link issue (use 'type=os.path.abspath' in -i argument). estrain@0: ############################begin the real analysis estrain@0: if analysis_mode=="a": estrain@0: if data_type in ["1","2","3"]:#use allele mode estrain@0: for_fq,rev_fq=get_input_files(make_dir,input_file,data_type,dirpath) estrain@0: os.chdir(make_dir) estrain@0: ###add a function to tell input files estrain@0: fnameA=for_fq.split("/")[-1] estrain@0: fnameB=rev_fq.split("/")[-1] estrain@0: current_time=time.strftime("%Y_%m_%d_%H_%M_%S", time.localtime()) estrain@0: sam,bam,sorted_bam,mapped_fq1,mapped_fq2,combined_fq,for_sai,rev_sai=get_temp_file_names(fnameA,fnameB) #get temp files id estrain@0: map_and_sort(threads,database,fnameA,fnameB,sam,bam,for_sai,rev_sai,sorted_bam,mapping_mode) #do mapping and sort estrain@0: xmlfile,new_fasta=extract_mapped_reads_and_do_assembly_and_blast(current_time,sorted_bam,combined_fq,mapped_fq1,mapped_fq2,threads,fnameA,fnameB,database,mapping_mode) #extract the mapped reads and do micro assembly and blast estrain@0: if xmlfile=="NA": estrain@0: O_choice,fliC_choice,fljB_choice,special_gene_list,contamination_O,contamination_H=("-","-","-",[],"","") estrain@0: else: estrain@0: Final_list=xml_parse_score_comparision_seqsero(xmlfile) #analyze xml and get parsed results estrain@0: file=open("data_log.txt","a") estrain@0: for x in Final_list: estrain@0: file.write("\t".join(str(y) for y in x)+"\n") estrain@0: file.close() estrain@0: Final_list_passed=[x for x in Final_list if float(x[0].split("_cov_")[1].split("_")[0])>=0.9 and (x[1]>=int(x[0].split("__")[1]) or x[1]>=int(x[0].split("___")[1].split("_")[3]) or x[1]>1000)] estrain@0: O_choice,fliC_choice,fljB_choice,special_gene_list,contamination_O,contamination_H,Otypes_uniq,H1_cont_stat_list,H2_cont_stat_list=predict_O_and_H_types(Final_list,Final_list_passed,new_fasta) #predict O, fliC and fljB estrain@0: subspecies=judge_subspecies(fnameA,dirpath) #predict subspecies estrain@0: ###output estrain@0: predict_form,predict_sero,star,star_line,claim=seqsero_from_formula_to_serotypes(O_choice,fliC_choice,fljB_choice,special_gene_list,subspecies) estrain@0: claim="" #04132019, disable claim for new report requirement estrain@0: contamination_report="" estrain@0: H_list=["fliC_"+x for x in H1_cont_stat_list if len(x)>0]+["fljB_"+x for x in H2_cont_stat_list if len(x)>0] estrain@0: if contamination_O!="" and contamination_H=="": estrain@0: contamination_report="#Potential inter-serotype contamination detected from O antigen signals. All O-antigens detected:"+"\t".join(Otypes_uniq)+"." estrain@0: elif contamination_O=="" and contamination_H!="": estrain@0: contamination_report="#Potential inter-serotype contamination detected or potential thrid H phase from H antigen signals. All H-antigens detected:"+"\t".join(H_list)+"." estrain@0: elif contamination_O!="" and contamination_H!="": estrain@0: contamination_report="#Potential inter-serotype contamination detected from both O and H antigen signals.All O-antigens detected:"+"\t".join(Otypes_uniq)+". All H-antigens detected:"+"\t".join(H_list)+"." estrain@0: if contamination_report!="": estrain@0: #contamination_report="potential inter-serotype contamination detected (please refer below antigen signal report for details)." #above contamination_reports are for back-up and bug fixing #web-based mode need to be re-used, 04132019 estrain@0: contamination_report=" Co-existence of multiple serotypes detected, indicating potential inter-serotype contamination. See 'Extracted_antigen_alleles.fasta' for detected serotype determinant alleles." estrain@0: #claim="\n"+open("Extracted_antigen_alleles.fasta","r").read()#used to store H and O antigen sequeences #04132019, need to change if using web-version estrain@0: ## ed_SL_09272019: change for new output format estrain@0: #if contamination_report+star_line+claim=="": #0413, new output style estrain@0: # note="" estrain@0: #else: estrain@0: # note="Note:" estrain@0: if clean_mode: estrain@0: subprocess.check_call("rm -rf ../"+make_dir,shell=True) estrain@0: make_dir="none-output-directory due to '-c' flag" estrain@0: else: estrain@0: new_file=open("SeqSero_result.txt","w") estrain@0: if O_choice=="": estrain@0: O_choice="-" estrain@0: if "N/A" not in predict_sero: estrain@0: new_file.write("Output_directory:\t"+make_dir+"\n"+ estrain@0: "Input files:\t"+"\t".join(input_file)+"\n"+ estrain@0: "O antigen prediction:\t"+O_choice+"\n"+ estrain@0: "H1 antigen prediction(fliC):\t"+fliC_choice+"\n"+ estrain@0: "H2 antigen prediction(fljB):\t"+fljB_choice+"\n"+ estrain@0: "Predicted subspecies:\t"+subspecies+"\n"+ estrain@0: "Predicted antigenic profile:\t"+predict_form+"\n"+ estrain@0: "Predicted serotype:\t"+predict_sero+"\n"+ # ed_SL_04152019: change serotype(s) to serotype estrain@0: note+contamination_report+star_line+claim+"\n")#+## estrain@0: else: estrain@0: #star_line=star_line.strip()+"\tNone such antigenic formula in KW.\n" estrain@0: star_line="" #04132019, for new output requirement, diable star_line if "NA" in output estrain@0: new_file.write("Output_directory:\t"+make_dir+"\n"+ estrain@0: "Input files:\t"+"\t".join(input_file)+"\n"+ estrain@0: "O antigen prediction:\t"+O_choice+"\n"+ estrain@0: "H1 antigen prediction(fliC):\t"+fliC_choice+"\n"+ estrain@0: "H2 antigen prediction(fljB):\t"+fljB_choice+"\n"+ estrain@0: "Predicted subspecies:\t"+subspecies+"\n"+ estrain@0: "Predicted antigenic profile:\t"+predict_form+"\n"+ estrain@0: "Predicted serotype:\t"+predict_form+"\n"+ # ed_SL_09242019: add serotype output for "N/A" prediction estrain@0: note+NA_note+contamination_report+star_line+claim+"\n")#+## estrain@0: new_file.close() estrain@0: print("\n") estrain@0: #subprocess.check_call("cat Seqsero_result.txt",shell=True) estrain@0: #subprocess.call("rm H_and_O_and_specific_genes.fasta* *.sra *.bam *.sam *.fastq *.gz *.fq temp.txt *.xml "+fnameA+"*_db* 2> /dev/null",shell=True) estrain@0: subprocess.call("rm H_and_O_and_specific_genes.fasta* *.sra *.bam *.sam *.fastq *.gz *.fq temp.txt "+fnameA+"*_db* 2> /dev/null",shell=True) estrain@0: if "N/A" not in predict_sero: estrain@0: #print("Output_directory:"+make_dir+"\nInput files:\t"+for_fq+" "+rev_fq+"\n"+"O antigen prediction:\t"+O_choice+"\n"+"H1 antigen prediction(fliC):\t"+fliC_choice+"\n"+"H2 antigen prediction(fljB):\t"+fljB_choice+"\n"+"Predicted antigenic profile:\t"+predict_form+"\n"+"Predicted subspecies:\t"+subspecies+"\n"+"Predicted serotype(s):\t"+predict_sero+star+"\nNote:"+contamination_report+star+star_line+claim+"\n")#+## estrain@0: print("Output_directory:\t"+make_dir+"\n"+ estrain@0: "Input files:\t"+"\t".join(input_file)+"\n"+ estrain@0: "O antigen prediction:\t"+O_choice+"\n"+ estrain@0: "H1 antigen prediction(fliC):\t"+fliC_choice+"\n"+ estrain@0: "H2 antigen prediction(fljB):\t"+fljB_choice+"\n"+ estrain@0: "Predicted subspecies:\t"+subspecies+"\n"+ estrain@0: "Predicted antigenic profile:\t"+predict_form+"\n"+ estrain@0: "Predicted serotype:\t"+predict_sero+"\n"+ # ed_SL_04152019: change serotype(s) to serotype estrain@0: note+contamination_report+star_line+claim+"\n")#+## estrain@0: else: estrain@0: print("Output_directory:\t"+make_dir+"\n"+ estrain@0: "Input files:\t"+"\t".join(input_file)+"\n"+ estrain@0: "O antigen prediction:\t"+O_choice+"\n"+ estrain@0: "H1 antigen prediction(fliC):\t"+fliC_choice+"\n"+ estrain@0: "H2 antigen prediction(fljB):\t"+fljB_choice+"\n"+ estrain@0: "Predicted subspecies:\t"+subspecies+"\n"+ estrain@0: "Predicted antigenic profile:\t"+predict_form+"\n"+ estrain@0: "Predicted serotype:\t"+predict_form+"\n"+ # ed_SL_09242019: add serotype output for "N/A" prediction estrain@0: note+NA_note+contamination_report+star_line+claim+"\n") estrain@0: else: estrain@0: print("Allele modes only support raw reads datatype, i.e. '-t 1 or 2 or 3'; please use '-m k'") estrain@0: elif analysis_mode=="k": estrain@0: #ex_dir = os.path.dirname(os.path.realpath(__file__)) estrain@0: ex_dir = os.path.abspath(os.path.join(os.path.dirname(os.path.dirname(__file__)),'seqsero2_db')) # ed_SL_09152019: change ex_dir for packaging estrain@0: #output_mode = args.mode estrain@0: for_fq,rev_fq=get_input_files(make_dir,input_file,data_type,dirpath) estrain@0: input_file = for_fq #-k will just use forward because not all reads were used estrain@0: os.chdir(make_dir) estrain@0: f = open(ex_dir + '/antigens.pickle', 'rb') estrain@0: lib_dict = pickle.load(f) estrain@0: f.close estrain@0: input_Ks=get_input_K(input_file,lib_dict,data_type,k_size) estrain@0: O_dict,H_dict,Special_dict=get_kmer_dict(lib_dict,input_Ks) estrain@0: highest_O,highest_fliC,highest_fljB=call_O_and_H_type(O_dict,H_dict,Special_dict,make_dir) estrain@0: subspecies=judge_subspecies_Kmer(Special_dict) estrain@0: if subspecies=="IIb" or subspecies=="IIa": estrain@0: subspecies="II" estrain@0: predict_form,predict_sero,star,star_line,claim = seqsero_from_formula_to_serotypes( estrain@0: highest_O.split('-')[1], highest_fliC, highest_fljB, Special_dict,subspecies) estrain@0: claim="" #no claim any more based on new output requirement estrain@0: ## ed_SL_09272019: change for new output format estrain@0: #if star_line+claim=="": #0413, new output style estrain@0: # note="" estrain@0: #else: estrain@0: # note="Note:" estrain@0: if clean_mode: estrain@0: subprocess.check_call("rm -rf ../"+make_dir,shell=True) estrain@0: make_dir="none-output-directory due to '-c' flag" estrain@0: ### ed_SL_05282019, fix the assignment issue of variable 'O_choice' using "-m k -c" estrain@0: if highest_O.split('-')[-1]=="": estrain@0: O_choice="-" estrain@0: else: estrain@0: O_choice=highest_O.split('-')[-1] estrain@0: ### estrain@0: else: estrain@0: if highest_O.split('-')[-1]=="": estrain@0: O_choice="-" estrain@0: else: estrain@0: O_choice=highest_O.split('-')[-1] estrain@0: #print("Output_directory:"+make_dir+"\tInput_file:"+input_file+"\tPredicted subpecies:"+subspecies + '\tPredicted antigenic profile:' + predict_form + '\tPredicted serotype(s):' + predict_sero) estrain@0: new_file=open("SeqSero_result.txt","w") estrain@0: #new_file.write("Output_directory:"+make_dir+"\nInput files:\t"+input_file+"\n"+"O antigen prediction:\t"+O_choice+"\n"+"H1 antigen prediction(fliC):\t"+highest_fliC+"\n"+"H2 antigen prediction(fljB):\t"+highest_fljB+"\n"+"Predicted antigenic profile:\t"+predict_form+"\n"+"Predicted subspecies:\t"+subspecies+"\n"+"Predicted serotype(s):\t"+predict_sero+star+"\n"+star+star_line+claim+"\n")#+## estrain@0: if "N/A" not in predict_sero: estrain@0: new_file.write("Output_directory:\t"+make_dir+"\n"+ estrain@0: "Input files:\t"+input_file+"\n"+ estrain@0: "O antigen prediction:\t"+O_choice+"\n"+ estrain@0: "H1 antigen prediction(fliC):\t"+highest_fliC+"\n"+ estrain@0: "H2 antigen prediction(fljB):\t"+highest_fljB+"\n"+ estrain@0: "Predicted subspecies:\t"+subspecies+"\n"+ estrain@0: "Predicted antigenic profile:\t"+predict_form+"\n"+ estrain@0: "Predicted serotype:\t"+predict_sero+"\n"+ # ed_SL_04152019: change serotype(s) to serotype estrain@0: note+star_line+claim+"\n")#+## estrain@0: else: estrain@0: #star_line=star_line.strip()+"\tNone such antigenic formula in KW.\n" estrain@0: star_line = "" #changed for new output requirement, 04132019 estrain@0: new_file.write("Output_directory:\t"+make_dir+"\n"+ estrain@0: "Input files:\t"+input_file+"\n"+ estrain@0: "O antigen prediction:\t"+O_choice+"\n"+ estrain@0: "H1 antigen prediction(fliC):\t"+highest_fliC+"\n"+ estrain@0: "H2 antigen prediction(fljB):\t"+highest_fljB+"\n"+ estrain@0: "Predicted subspecies:\t"+subspecies+"\n"+ estrain@0: "Predicted antigenic profile:\t"+predict_form+"\n"+ estrain@0: "Predicted serotype:\t"+predict_form+"\n"+ # ed_SL_09242019: add serotype output for "N/A" prediction estrain@0: note+NA_note+star_line+claim+"\n")#+## estrain@0: new_file.close() estrain@0: subprocess.call("rm *.fasta* *.fastq *.gz *.fq temp.txt *.sra 2> /dev/null",shell=True) estrain@0: if "N/A" not in predict_sero: estrain@0: print("Output_directory:\t"+make_dir+"\n"+ estrain@0: "Input files:\t"+input_file+"\n"+ estrain@0: "O antigen prediction:\t"+O_choice+"\n"+ estrain@0: "H1 antigen prediction(fliC):\t"+highest_fliC+"\n"+ estrain@0: "H2 antigen prediction(fljB):\t"+highest_fljB+"\n"+ estrain@0: "Predicted subspecies:\t"+subspecies+"\n"+ estrain@0: "Predicted antigenic profile:\t"+predict_form+"\n"+ estrain@0: "Predicted serotype:\t"+predict_sero+"\n"+ # ed_SL_04152019: change serotype(s) to serotype estrain@0: note+star_line+claim+"\n")#+## estrain@0: else: estrain@0: print("Output_directory:\t"+make_dir+"\n"+ estrain@0: "Input files:\t"+input_file+"\n"+ estrain@0: "O antigen prediction:\t"+O_choice+"\n"+ estrain@0: "H1 antigen prediction(fliC):\t"+highest_fliC+"\n"+ estrain@0: "H2 antigen prediction(fljB):\t"+highest_fljB+"\n"+ estrain@0: "Predicted subspecies:\t"+subspecies+"\n"+ estrain@0: "Predicted antigenic profile:\t"+predict_form+"\n"+ estrain@0: "Predicted serotype:\t"+predict_form+"\n"+ # ed_SL_09242019: add serotype output for "N/A" prediction estrain@0: note+NA_note+star_line+claim+"\n")#+## estrain@0: estrain@0: if __name__ == '__main__': estrain@0: main()