Louis BECQUEY

added cm coordinates and consensus 2D to align_colum

......@@ -68,12 +68,13 @@ To help you design your own SQL requests, we provide a description of the databa
* `puckering`: Conformation of the ribose cycle (10 classes depending on the phase_angle value)
## Table `align_column`, for positions in multiple sequence alignments
* `column_id`: A unique identifier
* `rfam_acc`: The family's MSA the column belongs to
* `index_ali`: Position of the column in the alignment (starts at 1)
* `cm_coord`: Position of the column in the Rfam covariance model of the family (starts at 1). The value is NULL in portions that are insertions compared to the model.
* `freq_A`, `freq_C`, `freq_G`, `freq_U`, `freq_other`: Nucleotide frequencies in the alignment at this position
* `gap_percent`: The frequencies of gaps at this position in the alignment (between 0.0 and 1.0)
* `consensus`: A consensus character (ACGUN or '-') summarizing the column, if we can. If >75% of the sequences are gaps at this position, the gap is picked as consensus. Otherwise, A/C/G/U is chosen if >50% of the non-gap positions are A/C/G/U. Otherwise, N is the consensus.
* `cons_sec_struct`: A consensus secondary structure for the RNAs of the family, obtained from the Infernal alignement. The structure is well-nested (no pseudoknots) and the possible symbols are '.' (unpaired) or '(' and ')' (paired). The field is NULL in portions that are insertions compared to the Rfam model of the family, meaning that their is no consensus on the structure.
There always is an entry, for each family (rfam_acc), with index_ali = 0; gap_percent = 1.0; and nucleotide frequencies set to 0.0. This entry is used when the nucleotide frequencies cannot be determined because of local alignment issues.
......
......@@ -17,12 +17,12 @@
* [DONE] Get filtered versions of the sequence alignments containing the 3D chains, publicly available for download
* [DONE] Get a consensus residue for each alignement column
* [DONE] Get an option to limit the number of cores
* [DONE] Move to SILVA LSU release 138.1
* [UPCOMING] Automated annotation of detected Recurrent Interaction Networks (RINs), see http://carnaval.lri.fr/ .
* [UPCOMING] Possibly, automated detection of HLs and ILs from the 3D Motif Atlas (BGSU). Maybe. Their own website already does the job.
* [UPCOMING] Weight sequences in alignment to give more importance to rarer sequences
* [UPCOMING] Give both gap_percent and insertion_gap_percent
* A field estimating the quality of the sequence alignment in table family.
* Possibly, more metrics about the alignments coming from Infernal.
* Move to SILVA LSU release 138
* Run cmscan ourselves from the NDB instead of using Rfam-PDB mappings ? (Iff this actually makes a real difference, untested yet)
* Use and save Infernal alignment bounds and truncation information
......
......@@ -862,10 +862,10 @@ class Downloader:
try:
print(f"Downloading {unit} from SILVA...", end='', flush=True)
if unit == "LSU":
subprocess.run(["wget", "-nv", "http://www.arb-silva.de/fileadmin/arb_web_db/release_132/ARB_files/SILVA_132_LSURef_07_12_17_opt.arb.gz",
subprocess.run(["wget", "-nv", "https://www.arb-silva.de/fileadmin/arb_web_db/release_138_1/ARB_files/SILVA_138.1_LSURef_opt.arb.gz",
"-O", path_to_seq_data + "realigned/LSU.arb.gz"])
else:
subprocess.run(["wget", "-nv", "http://www.arb-silva.de/fileadmin/silva_databases/release_138/ARB_files/SILVA_138_SSURef_05_01_20_opt.arb.gz",
subprocess.run(["wget", "-nv", "https://www.arb-silva.de/fileadmin/arb_web_db/release_138_1/ARB_files/SILVA_138.1_SSURef_opt.arb.gz",
"-O", path_to_seq_data + "realigned/SSU.arb.gz"])
except:
warn(f"Error downloading the {unit} database from SILVA", error=True)
......@@ -1379,11 +1379,11 @@ class Pipeline:
for r in results:
align = AlignIO.read(path_to_seq_data + "realigned/" + r[0] + "++.afa", "fasta")
nb_3d_chains = len([1 for r in align if '[' in r.id])
if r[0] in SSU_set: # SSU v138 is used
nb_homologs = 2225272 # source: https://www.arb-silva.de/documentation/release-138/
if r[0] in SSU_set: # SSU v138.1 is used
nb_homologs = 2224740 # source: https://www.arb-silva.de/documentation/release-1381/
nb_total_homol = nb_homologs + nb_3d_chains
elif r[0] in LSU_set: # LSU v132 is used
nb_homologs = 198843 # source: https://www.arb-silva.de/documentation/release-132/
elif r[0] in LSU_set: # LSU v138.1 is used
nb_homologs = 227331 # source: https://www.arb-silva.de/documentation/release-1381/
nb_total_homol = nb_homologs + nb_3d_chains
else:
nb_total_homol = len(align)
......@@ -1433,6 +1433,32 @@ class Pipeline:
p.join()
exit(1)
def extractCMs(self):
"""
Extracts Rfam CMs of the families for which we have 3D structures and compresses
them for later use with cmscan.
"""
# Retrieve the list of families and write them to a "keys file"
allfams = []
with sqlite3.connect(runDir+"/results/RNANet.db", timeout=10.0) as conn:
conn.execute('pragma journal_mode=wal')
allfams = sql_ask_database(conn, "SELECT DISTINCT rfam_acc FROM chain ORDER BY rfam_acc ASC;")
allfams = [ r[0]+"\n" for r in allfams ]
if not len(allfams):
return
with open(runDir + "/results/available_families.txt", "w") as f:
f.writelines(allfams)
# Extract the families from Rfam.cm
os.makedirs(runDir + "/results/cm/", exist_ok=True)
subprocess.run(["cmfetch", "-f", "-o", runDir + "/results/cm/RNANet.cm", path_to_seq_data + "Rfam.cm", runDir + "/results/available_families.txt"])
os.remove(runDir + "/results/available_families.txt")
# Compress the cm database for use with cmscan
subprocess.run(["rm", "-f", runDir + "/results/cm/RNANet.cm.i1p", runDir + "/results/cm/RNANet.cm.i1i", runDir + "/results/cm/RNANet.cm.i1m", runDir + "/results/cm/RNANet.cm.i1f"])
subprocess.run(["cmpress", runDir + "/results/cm/RNANet.cm"])
def output_results(self):
"""Produces CSV files, archive them, and additional metadata files
......@@ -1665,8 +1691,8 @@ def sql_define_tables(conn):
pair_count_tSS SMALLINT,
pair_count_other SMALLINT,
UNIQUE (structure_id, chain_name, rfam_acc),
FOREIGN KEY(structure_id) REFERENCES structure(pdb_id),
FOREIGN KEY(rfam_acc) REFERENCES family(rfam_acc)
FOREIGN KEY(structure_id) REFERENCES structure(pdb_id) ON DELETE CASCADE,
FOREIGN KEY(rfam_acc) REFERENCES family(rfam_acc) ON DELETE CASCADE
);
CREATE TABLE IF NOT EXISTS nucleotide (
chain_id INT,
......@@ -1721,6 +1747,7 @@ def sql_define_tables(conn):
CREATE TABLE IF NOT EXISTS align_column (
rfam_acc CHAR(7) NOT NULL,
index_ali INT NOT NULL,
cm_coord INT,
freq_A REAL,
freq_C REAL,
freq_G REAL,
......@@ -1728,12 +1755,17 @@ def sql_define_tables(conn):
freq_other REAL,
gap_percent REAL,
consensus CHAR(1),
cons_sec_struct CHAR(1),
PRIMARY KEY (rfam_acc, index_ali),
FOREIGN KEY(rfam_acc) REFERENCES family(rfam_acc)
FOREIGN KEY(rfam_acc) REFERENCES family(rfam_acc) ON DELETE CASCADE
);
""")
conn.commit()
# Prepare the WAL files while we're in single thread, otherwise it sometimes fails
# at the first access in WAL mode
conn.execute("pragma journal_mode=wal")
@trace_unhandled_exceptions
def sql_ask_database(conn, sql, warn_every=10):
"""
......@@ -1966,12 +1998,12 @@ def work_infer_mappings(update_only, allmappings, fullinference, codelist) -> li
sel_5_to_3 = (inferred_mappings['pdb_start'] < inferred_mappings['pdb_end'])
thisfam_5_3 = (inferred_mappings['rfam_acc'] == rfam) & sel_5_to_3
thisfam_3_5 = (inferred_mappings['rfam_acc'] == rfam) & (sel_5_to_3 == False)
print()
warn(f"Found mappings to {rfam} in both directions on the same interval, keeping only the 5'->3' one.")
# print()
# warn(f"Found mappings to {rfam} in both directions on the same interval, keeping only the 5'->3' one.")
else:
warn(f"There are mappings for {rfam} in both directions:", error=True)
print(inferred_mappings)
exit(1)
# exit(1)
# Compute consensus for chains in 5' -> 3' sense
if len(inferred_mappings[thisfam_5_3]):
......@@ -2297,8 +2329,10 @@ def work_realign(rfam_acc):
# Align the new sequences
with open(new_ali_path, 'w') as o:
p1 = subprocess.run(["cmalign", path_to_seq_data + f"realigned/{rfam_acc}.cm",
path_to_seq_data + f"realigned/{rfam_acc}_new.fa"],
p1 = subprocess.run(["cmalign", "--nonbanded", "--ifile", path_to_seq_data + f"realigned/{rfam_acc}.ins",
"--sfile", path_to_seq_data + f"realigned/{rfam_acc}.tsv",
path_to_seq_data + f"realigned/{rfam_acc}.cm",
path_to_seq_data + f"realigned/{rfam_acc}_new.fa"],
stdout=o, stderr=subprocess.PIPE)
notify("Aligned new sequences together")
......@@ -2337,8 +2371,8 @@ def work_realign(rfam_acc):
# Alignment does not exist yet. We need to compute it from scratch.
print(f"\t> Aligning {rfam_acc} sequences together (cmalign) ...", end='', flush=True)
p = subprocess.run(["cmalign", "--small", "--cyk", "--noprob", "--nonbanded", "--notrunc",
'-o', path_to_seq_data + f"realigned/{rfam_acc}++.stk",
p = subprocess.run(["cmalign", "--nonbanded", '-o', path_to_seq_data + f"realigned/{rfam_acc}++.stk",
"--ifile", path_to_seq_data + f"realigned/{rfam_acc}.ins", "--sfile", path_to_seq_data + f"realigned/{rfam_acc}.tsv",
path_to_seq_data + f"realigned/{rfam_acc}.cm",
path_to_seq_data + f"realigned/{rfam_acc}++.fa"],
stdout=subprocess.DEVNULL, stderr=subprocess.PIPE)
......@@ -2532,6 +2566,30 @@ def work_pssm_remap(f):
setproctitle(f"RNAnet.py work_pssm_remap({f}) saving")
# Get back the information of match/insertion states from the STK file
alignstk = AlignIO.read(path_to_seq_data + "realigned/" + f + "++.stk", "stockholm")
consensus_2d = alignstk.column_annotations["secondary_structure"]
del alignstk
cm_coord = 1
cm_coords = []
cm_2d = []
for x in consensus_2d:
if x == ".":
cm_coords.append(None)
cm_2d.append(None)
else:
cm_coords.append(cm_coord)
if x in "[(<{":
cm_2d.append("(")
elif x in "])>}":
cm_2d.append(")")
elif x in ",_-:":
cm_2d.append(".")
else:
warn("Unsupported WUSS secondary structure symbol : "+x)
cm_2d.append(".")
cm_coord += 1
# Save the re_mappings
conn = sqlite3.connect(runDir + '/results/RNANet.db', timeout=20.0)
conn.execute('pragma journal_mode=wal') # Allow multiple other readers to ask things while we execute this writing query
......@@ -2550,20 +2608,20 @@ def work_pssm_remap(f):
conn.commit()
# Save the useful columns in the database
data = [(f, j) + tuple(pssm_info[:,j-1]) + (consensus[j-1],) for j in sorted(columns_to_save)]
sql_execute(conn, """INSERT INTO align_column (rfam_acc, index_ali, freq_A, freq_C, freq_G, freq_U, freq_other, gap_percent, consensus)
VALUES (?, ?, ?, ?, ?, ?, ?, ?, ?) ON CONFLICT(rfam_acc, index_ali) DO
UPDATE SET freq_A=excluded.freq_A, freq_C=excluded.freq_C, freq_G=excluded.freq_G, freq_U=excluded.freq_U,
freq_other=excluded.freq_other, gap_percent=excluded.gap_percent, consensus=excluded.consensus;""", many=True, data=data)
data = [(f, j, cm_coords[j-1]) + tuple(pssm_info[:,j-1]) + (consensus[j-1], cm_2d[j-1]) for j in sorted(columns_to_save)]
sql_execute(conn, """INSERT INTO align_column (rfam_acc, index_ali, cm_coord, freq_A, freq_C, freq_G, freq_U, freq_other, gap_percent, consensus, cons_sec_struct)
VALUES (?, ?, ?, ?, ?, ?, ?, ?, ?, ?, ?) ON CONFLICT(rfam_acc, index_ali) DO
UPDATE SET cm_coord=excluded.cm_coord, freq_A=excluded.freq_A, freq_C=excluded.freq_C, freq_G=excluded.freq_G, freq_U=excluded.freq_U,
freq_other=excluded.freq_other, gap_percent=excluded.gap_percent, consensus=excluded.consensus, cons_sec_struct=excluded.cons_sec_struct;""", many=True, data=data)
# Add an unknown values column, with index_ali 0 (for nucleotides unsolved in 3D giving a gap '-' but found facing letter in the alignment)
sql_execute(conn, f"""INSERT OR IGNORE INTO align_column (rfam_acc, index_ali, freq_A, freq_C, freq_G, freq_U, freq_other, gap_percent, consensus)
VALUES (?, 0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, '-');""", data=(f,))
sql_execute(conn, f"""INSERT OR IGNORE INTO align_column (rfam_acc, index_ali, cm_coord, freq_A, freq_C, freq_G, freq_U, freq_other, gap_percent, consensus, cons_sec_struct)
VALUES (?, 0, NULL, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, '-', "NULL");""", data=(f,))
# Save the number of "used columns" to table family ( = the length of the alignment if it was composed only of the RNANet chains)
sql_execute(conn, f"UPDATE family SET ali_filtered_len = ? WHERE rfam_acc = ?;", data=(len(columns_to_save), f))
conn.close()
##########################################################################################
# Saving the filtered alignement with only the saved positinos
# Saving the filtered alignement with only the saved positions
##########################################################################################
setproctitle(f"RNAnet.py work_pssm_remap({f}) filtering alignment")
......@@ -2600,8 +2658,8 @@ def work_save(c, homology=True):
conn.execute('pragma journal_mode=wal')
if homology:
df = pd.read_sql_query(f"""
SELECT index_chain, old_nt_resnum, nt_position, nt_name, nt_code, nt_align_code,
is_A, is_C, is_G, is_U, is_other, freq_A, freq_C, freq_G, freq_U, freq_other, gap_percent, consensus, dbn,
SELECT index_chain, old_nt_resnum, nt_position, nt_name, nt_code, nt_align_code, cm_coord,
is_A, is_C, is_G, is_U, is_other, freq_A, freq_C, freq_G, freq_U, freq_other, gap_percent, consensus, cons_sec_struct, dbn,
paired, nb_interact, pair_type_LW, pair_type_DSSR, alpha, beta, gamma, delta, epsilon, zeta, epsilon_zeta,
chi, bb_type, glyco_bond, form, ssZp, Dp, eta, theta, eta_prime, theta_prime, eta_base, theta_base,
v0, v1, v2, v3, v4, amplitude, phase_angle, puckering FROM
......@@ -2682,6 +2740,7 @@ if __name__ == "__main__":
work_save(c, homology=False)
print("Completed.")
exit(0)
# At this point, structure, chain and nucleotide tables of the database are up to date.
# (Modulo some statistics computed by statistics.py)
......@@ -2716,6 +2775,7 @@ if __name__ == "__main__":
idxQueue = thr_idx_mgr.Queue()
pp.remap()
pp.extractCMs()
# At this point, the align_column and re_mapping tables are up-to-date.
......
2z9q_1_A_1-72
1ml5_1_b_5-121
1ml5_1_a_1-2914
3ep2_1_Y_1-72
3eq3_1_Y_1-72
4v48_1_A6_1-73
1ml5_1_A_2-1520
1qzb_1_B_1-73
1qza_1_B_1-73
1ls2_1_B_1-73
1gsg_1_T_1-72
3jcr_1_H_1-115
1vy7_1_AY_1-73
1vy7_1_CY_1-73
4w2h_1_CY_1-73
5zzm_1_M_3-118
2rdo_1_A_3-118
4v48_1_A9_3-118
4v47_1_A9_3-118
2ob7_1_A_10-319
1x1l_1_A_1-132
1zc8_1_Z_1-93
2ob7_1_D_1-132
4v42_1_BB_5-121
4v42_1_BA_1-2914
1r2x_1_C_1-58
1r2w_1_C_1-58
1eg0_1_L_1-56
5zzm_1_N_1-2904
2rdo_1_B_1-2904
3dg2_1_B_1-2904
3dg0_1_B_1-2904
4v48_1_A0_1-2904
4v47_1_A0_1-2904
3dg4_1_B_1-2904
3dg5_1_B_1-2904
3dg2_1_A_1-1542
3dg0_1_A_1-1542
4v48_1_BA_1-1543
4v47_1_BA_1-1542
3dg4_1_A_1-1542
3dg5_1_A_1-1542
1eg0_1_O_1-73
1zc8_1_A_1-59
1mvr_1_D_1-61
4adx_1_9_1-123
1zn1_1_B_1-59
1jgq_1_A_2-1520
4v42_1_AA_2-1520
1jgo_1_A_2-1520
1jgp_1_A_2-1520
1emi_1_B_1-108
3iy9_1_A_498-1027
3ep2_1_B_1-50
3eq3_1_B_1-50
3eq4_1_B_1-50
3pgw_1_R_1-164
3pgw_1_N_1-164
3cw1_1_x_1-138
3cw1_1_w_1-138
3cw1_1_V_1-138
3cw1_1_v_1-138
2iy3_1_B_9-105
3jcr_1_N_1-107
2vaz_1_A_64-177
2ftc_1_R_81-1466
3jcr_1_M_1-141
4v5z_1_B0_1-2902
5g2x_1_A_595-692
3iy8_1_A_1-540
4v5z_1_BY_2-113
4v5z_1_BZ_1-70
4v5z_1_B1_2-125
4adx_1_0_1-2925
1mvr_1_B_3-96
3eq4_1_Y_1-69
6uz7_1_8_2140-2827
4v5z_1_AA_1-1563
2z9q_1_A_1-72
DSSR warning 2z9q.json: no nucleotides found. Ignoring 2z9q_1_A_1-72.
1ml5_1_b_5-121
Could not find nucleotides of chain b in annotation 1ml5.json. Either there is a problem with 1ml5 mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
1ml5_1_a_1-2914
Could not find nucleotides of chain a in annotation 1ml5.json. Either there is a problem with 1ml5 mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
3ep2_1_Y_1-72
DSSR warning 3ep2.json: no nucleotides found. Ignoring 3ep2_1_Y_1-72.
3eq3_1_Y_1-72
DSSR warning 3eq3.json: no nucleotides found. Ignoring 3eq3_1_Y_1-72.
4v48_1_A6_1-73
DSSR warning 4v48.json: no nucleotides found. Ignoring 4v48_1_A6_1-73.
1ml5_1_A_2-1520
Could not find nucleotides of chain A in annotation 1ml5.json. Either there is a problem with 1ml5 mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
1qzb_1_B_1-73
DSSR warning 1qzb.json: no nucleotides found. Ignoring 1qzb_1_B_1-73.
1qza_1_B_1-73
DSSR warning 1qza.json: no nucleotides found. Ignoring 1qza_1_B_1-73.
1ls2_1_B_1-73
DSSR warning 1ls2.json: no nucleotides found. Ignoring 1ls2_1_B_1-73.
1gsg_1_T_1-72
DSSR warning 1gsg.json: no nucleotides found. Ignoring 1gsg_1_T_1-72.
3jcr_1_H_1-115
DSSR warning 3jcr.json: no nucleotides found. Ignoring 3jcr_1_H_1-115.
1vy7_1_AY_1-73
Sequence is too short. (< 5 resolved nts)
1vy7_1_CY_1-73
Sequence is too short. (< 5 resolved nts)
4w2h_1_CY_1-73
Sequence is too short. (< 5 resolved nts)
5zzm_1_M_3-118
DSSR warning 5zzm.json: no nucleotides found. Ignoring 5zzm_1_M_3-118.
2rdo_1_A_3-118
DSSR warning 2rdo.json: no nucleotides found. Ignoring 2rdo_1_A_3-118.
4v48_1_A9_3-118
DSSR warning 4v48.json: no nucleotides found. Ignoring 4v48_1_A9_3-118.
4v47_1_A9_3-118
DSSR warning 4v47.json: no nucleotides found. Ignoring 4v47_1_A9_3-118.
2ob7_1_A_10-319
DSSR warning 2ob7.json: no nucleotides found. Ignoring 2ob7_1_A_10-319.
1x1l_1_A_1-132
DSSR warning 1x1l.json: no nucleotides found. Ignoring 1x1l_1_A_1-132.
1zc8_1_Z_1-93
DSSR warning 1zc8.json: no nucleotides found. Ignoring 1zc8_1_Z_1-93.
2ob7_1_D_1-132
DSSR warning 2ob7.json: no nucleotides found. Ignoring 2ob7_1_D_1-132.
4v42_1_BB_5-121
Could not find nucleotides of chain BB in annotation 4v42.json. Either there is a problem with 4v42 mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
4v42_1_BA_1-2914
Could not find nucleotides of chain BA in annotation 4v42.json. Either there is a problem with 4v42 mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
1r2x_1_C_1-58
DSSR warning 1r2x.json: no nucleotides found. Ignoring 1r2x_1_C_1-58.
1r2w_1_C_1-58
DSSR warning 1r2w.json: no nucleotides found. Ignoring 1r2w_1_C_1-58.
1eg0_1_L_1-56
DSSR warning 1eg0.json: no nucleotides found. Ignoring 1eg0_1_L_1-56.
5zzm_1_N_1-2904
DSSR warning 5zzm.json: no nucleotides found. Ignoring 5zzm_1_N_1-2904.
2rdo_1_B_1-2904
DSSR warning 2rdo.json: no nucleotides found. Ignoring 2rdo_1_B_1-2904.
3dg2_1_B_1-2904
DSSR warning 3dg2.json: no nucleotides found. Ignoring 3dg2_1_B_1-2904.
3dg0_1_B_1-2904
DSSR warning 3dg0.json: no nucleotides found. Ignoring 3dg0_1_B_1-2904.
4v48_1_A0_1-2904
DSSR warning 4v48.json: no nucleotides found. Ignoring 4v48_1_A0_1-2904.
4v47_1_A0_1-2904
DSSR warning 4v47.json: no nucleotides found. Ignoring 4v47_1_A0_1-2904.
3dg4_1_B_1-2904
DSSR warning 3dg4.json: no nucleotides found. Ignoring 3dg4_1_B_1-2904.
3dg5_1_B_1-2904
DSSR warning 3dg5.json: no nucleotides found. Ignoring 3dg5_1_B_1-2904.
3dg2_1_A_1-1542
DSSR warning 3dg2.json: no nucleotides found. Ignoring 3dg2_1_A_1-1542.
3dg0_1_A_1-1542
DSSR warning 3dg0.json: no nucleotides found. Ignoring 3dg0_1_A_1-1542.
4v48_1_BA_1-1543
DSSR warning 4v48.json: no nucleotides found. Ignoring 4v48_1_BA_1-1543.
4v47_1_BA_1-1542
DSSR warning 4v47.json: no nucleotides found. Ignoring 4v47_1_BA_1-1542.
3dg4_1_A_1-1542
DSSR warning 3dg4.json: no nucleotides found. Ignoring 3dg4_1_A_1-1542.
3dg5_1_A_1-1542
DSSR warning 3dg5.json: no nucleotides found. Ignoring 3dg5_1_A_1-1542.
1eg0_1_O_1-73
DSSR warning 1eg0.json: no nucleotides found. Ignoring 1eg0_1_O_1-73.
1zc8_1_A_1-59
DSSR warning 1zc8.json: no nucleotides found. Ignoring 1zc8_1_A_1-59.
1mvr_1_D_1-61
DSSR warning 1mvr.json: no nucleotides found. Ignoring 1mvr_1_D_1-61.
4adx_1_9_1-123
DSSR warning 4adx.json: no nucleotides found. Ignoring 4adx_1_9_1-123.
1zn1_1_B_1-59
DSSR warning 1zn1.json: no nucleotides found. Ignoring 1zn1_1_B_1-59.
1jgq_1_A_2-1520
Could not find nucleotides of chain A in annotation 1jgq.json. Either there is a problem with 1jgq mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
4v42_1_AA_2-1520
Could not find nucleotides of chain AA in annotation 4v42.json. Either there is a problem with 4v42 mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
1jgo_1_A_2-1520
Could not find nucleotides of chain A in annotation 1jgo.json. Either there is a problem with 1jgo mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
1jgp_1_A_2-1520
Could not find nucleotides of chain A in annotation 1jgp.json. Either there is a problem with 1jgp mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
1emi_1_B_1-108
DSSR warning 1emi.json: no nucleotides found. Ignoring 1emi_1_B_1-108.
3iy9_1_A_498-1027
DSSR warning 3iy9.json: no nucleotides found. Ignoring 3iy9_1_A_498-1027.
3ep2_1_B_1-50
DSSR warning 3ep2.json: no nucleotides found. Ignoring 3ep2_1_B_1-50.
3eq3_1_B_1-50
DSSR warning 3eq3.json: no nucleotides found. Ignoring 3eq3_1_B_1-50.
3eq4_1_B_1-50
DSSR warning 3eq4.json: no nucleotides found. Ignoring 3eq4_1_B_1-50.
3pgw_1_R_1-164
DSSR warning 3pgw.json: no nucleotides found. Ignoring 3pgw_1_R_1-164.
3pgw_1_N_1-164
DSSR warning 3pgw.json: no nucleotides found. Ignoring 3pgw_1_N_1-164.
3cw1_1_x_1-138
DSSR warning 3cw1.json: no nucleotides found. Ignoring 3cw1_1_x_1-138.
3cw1_1_w_1-138
DSSR warning 3cw1.json: no nucleotides found. Ignoring 3cw1_1_w_1-138.
3cw1_1_V_1-138
DSSR warning 3cw1.json: no nucleotides found. Ignoring 3cw1_1_V_1-138.
3cw1_1_v_1-138
DSSR warning 3cw1.json: no nucleotides found. Ignoring 3cw1_1_v_1-138.
2iy3_1_B_9-105
DSSR warning 2iy3.json: no nucleotides found. Ignoring 2iy3_1_B_9-105.
3jcr_1_N_1-107
DSSR warning 3jcr.json: no nucleotides found. Ignoring 3jcr_1_N_1-107.
2vaz_1_A_64-177
DSSR warning 2vaz.json: no nucleotides found. Ignoring 2vaz_1_A_64-177.
2ftc_1_R_81-1466
DSSR warning 2ftc.json: no nucleotides found. Ignoring 2ftc_1_R_81-1466.
3jcr_1_M_1-141
DSSR warning 3jcr.json: no nucleotides found. Ignoring 3jcr_1_M_1-141.
4v5z_1_B0_1-2902
DSSR warning 4v5z.json: no nucleotides found. Ignoring 4v5z_1_B0_1-2902.
5g2x_1_A_595-692
Sequence is too short. (< 5 resolved nts)
3iy8_1_A_1-540
DSSR warning 3iy8.json: no nucleotides found. Ignoring 3iy8_1_A_1-540.
4v5z_1_BY_2-113
DSSR warning 4v5z.json: no nucleotides found. Ignoring 4v5z_1_BY_2-113.
4v5z_1_BZ_1-70
DSSR warning 4v5z.json: no nucleotides found. Ignoring 4v5z_1_BZ_1-70.
4v5z_1_B1_2-125
DSSR warning 4v5z.json: no nucleotides found. Ignoring 4v5z_1_B1_2-125.
4adx_1_0_1-2925
DSSR warning 4adx.json: no nucleotides found. Ignoring 4adx_1_0_1-2925.
1mvr_1_B_3-96
DSSR warning 1mvr.json: no nucleotides found. Ignoring 1mvr_1_B_3-96.
3eq4_1_Y_1-69
DSSR warning 3eq4.json: no nucleotides found. Ignoring 3eq4_1_Y_1-69.
6uz7_1_8_2140-2827
Could not find nucleotides of chain 8 in annotation 6uz7.json. Either there is a problem with 6uz7 mmCIF download, or the bases are not resolved in the structure. Delete it and retry.
4v5z_1_AA_1-1563
DSSR warning 4v5z.json: no nucleotides found. Ignoring 4v5z_1_AA_1-1563.
......@@ -958,6 +958,8 @@ def par_distance_matrix(filelist, f, label, consider_all_atoms, s):
@trace_unhandled_exceptions
def get_avg_std_distance_matrix(f, consider_all_atoms, multithread=False):
np.seterr(divide='ignore') # ignore division by zero issues
if consider_all_atoms:
label = "base"
else:
......@@ -1024,44 +1026,38 @@ def get_avg_std_distance_matrix(f, consider_all_atoms, multithread=False):
p.join()
exit(1)
with warnings.catch_warnings(record=True) as w:
warnings.simplefilter("error")
# Calculation of the average matrix
try:
avg = avg/counts
assert (counts >0).all()
except RuntimeWarning as e:
# there will be division by 0 if count is 0 at some position = gap only column
pass
np.savetxt(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '_average.csv' , avg, delimiter=",", fmt="%.3f")
fig, ax = plt.subplots()
im = ax.imshow(avg)
cbar = ax.figure.colorbar(im, ax=ax)
cbar.ax.set_ylabel("Angströms", rotation=-90, va="bottom")
ax.set_title(f"Average distance between {f} residues (Angströms)")
fig.tight_layout()
fig.savefig(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '_average.png', dpi=300)
plt.close()
# Calculation of the average matrix
avg = np.divide(avg, counts, where=counts>0, out=np.full_like(avg, np.NaN)) # Ultrafancy way to take avg/counts or NaN if counts is 0
np.savetxt(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '_average.csv' , avg, delimiter=",", fmt="%.3f")
fig, ax = plt.subplots()
im = ax.imshow(avg)
cbar = ax.figure.colorbar(im, ax=ax)
cbar.ax.set_ylabel("Angströms", rotation=-90, va="bottom")
ax.set_title(f"Average distance between {f} residues (Angströms)")
fig.tight_layout()
fig.savefig(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '_average.png', dpi=300)
plt.close()
# Calculation of the standard deviation matrix by the Huygens theorem
try:
std = np.sqrt(std/counts - np.power(avg/counts, 2))
assert (counts >0).all()
except RuntimeWarning as e:
# there will be division by 0 if count is 0 at some position = gap only column
pass
np.savetxt(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '_stdev.csv' , std, delimiter=",", fmt="%.3f")
fig, ax = plt.subplots()
im = ax.imshow(std)
cbar = ax.figure.colorbar(im, ax=ax)
cbar.ax.set_ylabel("Angströms", rotation=-90, va="bottom")
ax.set_title(f"Standard deviation of distances between {f} residues (Angströms)")
fig.tight_layout()
fig.savefig(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '_std.png', dpi=300)
plt.close()
# Calculation of the standard deviation matrix by the Huygens theorem
std = np.divide(std, counts, where=counts>0, out=np.full_like(std, np.NaN))
mask = np.invert(np.isnan(std))
value = std[mask] - np.power(avg[mask], 2)
if ((value[value<0] < -1e-2).any()):
warn("Erasing very negative variance value !")
value[value<0] = 0.0 # floating point problems !
std[mask] = np.sqrt(value)
np.savetxt(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '_stdev.csv' , std, delimiter=",", fmt="%.3f")
fig, ax = plt.subplots()
im = ax.imshow(std)
cbar = ax.figure.colorbar(im, ax=ax)
cbar.ax.set_ylabel("Angströms", rotation=-90, va="bottom")
ax.set_title(f"Standard deviation of distances between {f} residues (Angströms)")
fig.tight_layout()
fig.savefig(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '_std.png', dpi=300)
plt.close()
# Save log
with open(runDir + '/results/distance_matrices/' + f + '_'+ label + '/' + f + '.log', 'a') as logfile:
......