PDB search¶
This interactive Python notebook will guide you through various ways of programmatically accessing Protein Data Bank in Europe (PDBe) data using REST API
The REST API is a programmatic way to obtain information from the PDB and EMDB. You can access details about:
sample
experiment
models
compounds
cross-references
publications
quality
assemblies and more… For more information, visit https://www.ebi.ac.uk/pdbe/pdbe-rest-api
This notebook is a part of the training material series, and focuses on getting information from the PDBe search API. Retrieve this material and many more from GitHub
1) Making imports and setting variables¶
First, we import some packages that we will use, and set some variables.
Note: Full list of valid URLs is available from https://www.ebi.ac.uk/pdbe/api/doc/
[162]:
import requests # used for getting data from a URL
from pprint import pprint # pretty print
import pandas as pd # used for turning results into mini databases
from solrq import Q # used to turn result queries into the right format
import cufflinks as cf
import plotly.offline as py
cf.go_offline() # required to use plotly offline (no account required).
py.init_notebook_mode() # graphs charts inline (IPython).
search_url = "https://www.ebi.ac.uk/pdbe/search/pdb/select?" # the rest of the URL used for PDBe's search API.
2) a function to get data from the search API¶
Let’s start with defining a function that can be used to GET data from the PDBe search API.
[163]:
def make_request(search_dict, number_of_rows=10):
"""
makes a get request to the PDBe API
:param dict search_dict: the terms used to search
:param number_of_rows: number or rows to return - limited to 10
:return dict: response JSON
"""
if 'rows' not in search_dict:
search_dict['rows'] = number_of_rows
search_dict['wt'] = 'json'
# pprint(search_dict)
response = requests.post(search_url, data=search_dict)
if response.status_code == 200:
return response.json()
else:
print("[No data retrieved - %s] %s" % (response.status_code, response.text))
return {}
3) formatting the search terms¶
This will allow us to use human readable search terms and this function will make a URL that the search API can handle.
[164]:
def format_search_terms(search_terms, filter_terms=None):
ret = {'q': str(search_terms)}
if filter_terms:
fl = '{}'.format(','.join(filter_terms))
ret['fl'] = fl
return ret
4) Getting useful data out of the search¶
This function will run the search and will return a list of the results
[165]:
def run_search(search_terms, filter_terms=None, number_of_rows=100):
search_term = format_search_terms(search_terms, filter_terms)
response = make_request(search_term, number_of_rows)
results = response.get('response', {}).get('docs', [])
print('Number of results: {}'.format(len(results)))
return results
5) running a search¶
Now we are ready to actually run a search against the PDB API for entries containing human Dihydrofolate reductase in the PDB. This will return a list of results - only 10 to start with.
A list of search terms is available at: https://www.ebi.ac.uk/pdbe/api/doc/search.html
This will return details of human Dihydrofolate reductase’s in the PDB
We are going to use the Python package “solrq” to make sure we get the search in the right format.
Q(molecule_name=“Dihydrofolate reductase”)
Here we are searching for molecules named Dihydrofolate reductase. If we search for two terms i.e. molecule_name and organism_scientific_name then we will get molecules that match both search terms.
We will return the number of results for two searches.
The first one will hit the limit of 100. There are more than 100 Dihydrofolate reductase structures. We have to add the argument “number_of_rows” to a higher number, say 1000, to find all the examples.
[166]:
print('1st search')
search_terms = Q(molecule_name="Dihydrofolate reductase")
first_results = run_search(search_terms)
1st search
Number of results: 100
[167]:
print('1st search - more rows')
search_terms = Q(molecule_name="Dihydrofolate reductase")
first_results_more_rows = run_search(search_terms, number_of_rows=1000)
1st search - more rows
Number of results: 405
We will add organism_name of Human to the query to limit the results to only return those that are structures of Human Dihydrofolate reductase.
[168]:
print('2nd search')
search_terms = Q(molecule_name="Dihydrofolate reductase",organism_name="Human")
second_results = run_search(search_terms)
2nd search
Number of results: 79
How did we know which search terms to use?
We will then look at the last result. We will print the data we have for the first result.
This will be the first item of the list “second_results” i.e. second_results[0]
We are using “pprint” (pretty print) rather than “print” to make the result easier to read.
All of the “keys” on the left side of the results can be used as a search term.
[169]:
pprint(second_results[0])
L',
'q_struct_asym_id': ['A'],
'q_structure_determination_method': ['MOLECULAR REPLACEMENT'],
'q_structure_solution_software': ['MOLREP'],
'q_superkingdom': ['Eukaryota'],
'q_tax_id': [9606],
'q_tax_query': [9606],
'q_title': 'Preferential Selection of Isomer Binding from Chiral Mixtures: '
'Alternate Binding Modes Observed for the E- and Z-isomers of a '
'Series of 5-Substituted 2,4-Diaminofuro[2,3-d]pyrimidines as '
'Ternary Complexes with NADPH and Human Dihydrofolate Reductase',
'q_uniprot': ['P00374 : DYR_HUMAN'],
'q_uniprot_accession': ['P00374', 'P00374-2'],
'q_uniprot_accession_best': ['P00374'],
'q_uniprot_best': ['P00374 : DYR_HUMAN'],
'q_uniprot_coverage': [0.99],
'q_uniprot_features': ['Protein has possible alternate isoforms',
'DHFR',
'Protein has possible natural variant ',
'Nucleotide binding - NADP'],
'q_uniprot_id': ['DYR_HUMAN', 'DYR_HUMAN'],
'q_uniprot_id_best': ['DYR_HUMAN'],
'q_uniprot_non_canonical': ['P00374-2 : DYR_HUMAN'],
'q_unp_count': 1,
'q_unp_nf90_accession': ['B0YJ76',
'A0A2K6C3Y8',
'S5WD14',
'S5VM81',
'P00374',
'A0A024RAQ3'],
'q_unp_nf90_id': ['B0YJ76_HUMAN',
'A0A2K6C3Y8_MACNE',
'S5WD14_SHISS',
'S5VM81_ECO57',
'DYR_HUMAN',
'A0A024RAQ3_HUMAN'],
'q_unp_nf90_organism': ['Homo sapiens (Human)',
'Macaca nemestrina (Pig-tailed macaque)',
'Shigella sonnei (strain Ss046)',
'Escherichia coli O157:H7',
'Homo sapiens (Human)',
'Homo sapiens (Human)'],
'q_unp_nf90_protein_name': ['Dihydrofolate reductase',
'DHFR domain-containing protein',
'Dihydrofolate reductase',
'Dihydrofolate reductase',
'Dihydrofolate reductase',
'Dihydrofolate reductase'],
'q_unp_nf90_tax_id': ['9606', '9545', '300269', '83334', '9606', '9606'],
'r_factor': 0.19536,
'r_free': 0.2456,
'r_work': [0.19278],
'rank': ['species',
'genus',
'subfamily',
'family',
'order',
'class',
'phylum',
'kingdom',
'superkingdom',
'species',
'genus',
'family',
'order',
'class',
'phylum',
'superkingdom'],
'reactant_id': ['NDP'],
'refinement_software': ['REFMAC'],
'release_date': '2010-12-01T01:00:00Z',
'release_year': 2010,
'resolution': 1.9,
'revision_date': '2017-11-08T01:00:00Z',
'revision_year': 2017,
'sample_preparation_method': ['Engineered'],
'seq_100_cluster_number': '24820',
'seq_100_cluster_rank': 39,
'seq_30_cluster_number': '162',
'seq_30_cluster_rank': 65,
'seq_40_cluster_number': '16985',
'seq_40_cluster_rank': 65,
'seq_50_cluster_number': '28892',
'seq_50_cluster_rank': 65,
'seq_70_cluster_number': '25856',
'seq_70_cluster_rank': 65,
'seq_90_cluster_number': '15408',
'seq_90_cluster_rank': 57,
'seq_95_cluster_number': '4071',
'seq_95_cluster_rank': 51,
'serial_crystal_experiment': ['N'],
'spacegroup': 'H 3',
'status': 'REL',
'struct_asym_id': ['A'],
'structure_determination_method': ['MOLECULAR REPLACEMENT'],
'structure_solution_software': ['MOLREP'],
'superkingdom': ['Eukaryota'],
't_abstracttext_unassigned': ['The crystal structures of six human '
'dihydrofolate reductase (hDHFR) ternary '
'complexes with NADPH and a series of mixed E/Z '
'isomers of 5-substituted '
'5-[2-(2-methoxyphenyl)-prop-1-en-1-yl]furo[2,3-d]pyrimidine-2,4-diamines '
'substituted at the C9 position with propyl, '
'isopropyl, cyclopropyl, butyl, isobutyl and '
'sec-butyl (E2-E7, Z3) were determined and the '
'results were compared with the resolved E and '
'Z isomers of the C9-methyl parent compound. '
'The configuration of all of the inhibitors, '
'save one, was observed as the E isomer, in '
'which the binding of the furopyrimidine ring '
'is flipped such that the 4-amino group binds '
'in the 4-oxo site of folate. The Z3\xa0isomer '
'of the C9-isopropyl analog has the normal '
'2,4-diaminopyrimidine ring binding geometry, '
'with the furo oxygen near Glu30 and the '
'4-amino group interacting near the cofactor '
'nicotinamide ring. Electron-density maps for '
'these structures revealed the binding of only '
'one isomer to hDHFR, despite the fact that '
'chiral mixtures (E:Z ratios of 2:1, 3:1 and '
'3:2) of the inhibitors were incubated with '
'hDHFR prior to crystallization. Superposition '
'of the hDHFR complexes with E2 and Z3 shows '
"that the 2'-methoxyphenyl ring of E2 is "
'perpendicular to that of Z3. The most potent '
'inhibitor in this series is the isopropyl '
'analog Z3 and the least potent is the isobutyl '
'analog E6, consistent with data that show that '
'the Z isomer makes the most favorable '
'interactions with the active-site residues. '
'The isobutyl moiety of E6 is observed in two '
'orientations and the resultant steric crowding '
'of the E6 analog is consistent with its weaker '
'activity. The alternative binding modes '
'observed for the furopyrimidine ring in these '
'E/Z isomers suggest that new templates can be '
'designed to probe these binding regions of the '
'DHFR active site.'],
't_all_compound_names': ['Nicotinamide-adenine dinucleotide',
'NDP',
'D2F : '
'5-[(1E)-2-(2-methoxyphenyl)hex-1-en-1-yl]furo[2,3-d]pyrimidine-2,4-diamine',
'NDP : NADPH '
'DIHYDRO-NICOTINAMIDE-ADENINE-DINUCLEOTIDE PHOSPHATE',
'SO4 : SULFATE ION',
'NDP : Dihydronicotinamide-adenine dinucleotide '
'phosphate',
'NDP : Reduced nicotinamide adenine dinucleotide '
'phosphate',
'NDP : TPNH',
'SO4 : Sulfate',
'SO4 : Sulfate dianion',
'SO4 : Sulfate(2-)',
'SO4 : Sulfuric acid ion(2-)',
'D2F : '
'5-[(1E)-2-(2-methoxyphenyl)hex-1-en-1-yl]furo[2,3-d]pyrimidine-2,4-diamine',
'SO4 : sulfate',
'D2F : '
'5-[(E)-2-(2-methoxyphenyl)hex-1-enyl]furo[2,3-d]pyrimidine-2,4-diamine',
'NDP : '
'[[(2R,3S,4R,5R)-5-(3-aminocarbonyl-4H-pyridin-1-yl)-3,4-dihydroxy-oxolan-2-yl]methoxy-hydroxy-phosphoryl] '
'[(2R,3R,4R,5R)-5-(6-aminopurin-9-yl)-3-hydroxy-4-phosphonooxy-oxolan-2-yl]methyl '
'hydrogen phosphate'],
't_all_enzyme_names': ['Oxidoreductases',
'Acting on the CH-NH group of donors',
'With NAD(+) or NADP(+) as acceptor',
'Dihydrofolate reductase',
'1.5.1.3 : Dihydrofolate reductase',
'5,6,7,8-tetrahydrofolate:NADP(+) oxidoreductase'],
't_all_go_terms': ['cytoplasm',
'mitochondrion',
'cytosol',
'folic acid binding',
'oxidoreductase activity',
'NADPH binding',
'RNA binding',
'sequence-specific mRNA binding',
'mRNA binding',
'dihydrofolate reductase activity',
'NADP binding',
'methotrexate binding',
'translation repressor activity, mRNA regulatory element '
'binding',
'drug binding',
'one-carbon metabolic process',
'negative regulation of translation',
'folic acid metabolic process',
'response to methotrexate',
'regulation of removal of superoxide radicals',
'tetrahydrobiopterin biosynthetic process',
'tetrahydrofolate metabolic process',
'tetrahydrofolate biosynthetic process',
'oxidation-reduction process',
'regulation of transcription involved in G1/S transition '
'of mitotic cell cycle',
'positive regulation of nitric-oxide synthase activity',
'axon regeneration',
'dihydrofolate metabolic process'],
't_all_sequence_family': ['IPR001796 : Dihydrofolate reductase domain',
'IPR024072 : Dihydrofolate reductase-like domain '
'superfamily',
'PF00186 : DHFR_1',
'CL0387 : DHFred'],
't_all_structure_family': ['3-Layer(aba) Sandwich',
'Alpha Beta',
'3.40.430.10',
'Dihydrofolate Reductase, subunit A',
'Dihydrofolate Reductase, subunit A'],
't_citation_authors': ['Cody V', 'Piraino J', 'Pace J', 'Li W', 'Gangjee A'],
't_citation_title': ['Preferential selection of isomer binding from chiral '
'mixtures: alternate binding modes observed for the E '
'and Z isomers of a series of 5-substituted '
'2,4-diaminofuro[2,3-d]pyrimidines as ternary complexes '
'with NADPH and human dihydrofolate reductase.'],
't_entry_authors': ['Cody V'],
't_entry_info': ['POTASSIUM PHOSPHATE',
'AMMONIUM SULFATE',
'ETHANOL',
'MOSFLM',
'SCALA',
'Image plate',
'RIGAKU RAXIS IV',
'X-ray diffraction',
'REFMAC',
'MOLECULAR REPLACEMENT',
'MOLREP'],
't_entry_title': ['Preferential Selection of Isomer Binding from Chiral '
'Mixtures: Alternate Binding Modes Observed for the E- and '
'Z-isomers of a Series of 5-Substituted '
'2,4-Diaminofuro[2,3-d]pyrimidines as Ternary Complexes '
'with NADPH and Human Dihydrofolate Reductase'],
't_expression_organism_name': ['Escherichia coli',
'Bacterium Coli',
'Enterococcus Coli',
'Bacterium 10a',
'Escherichia Coli',
'Escherichia Sp. 3_2_53faa',
'Escherichia/Shigella Coli',
'Ecolx',
'Bacterium E3',
'Bacterium Coli Commune',
'E. Coli',
'Bacillus Coli',
'Escherichia Sp. Mar',
'Escherichia coli',
'Escherichia',
'Enterobacteriaceae',
'Enterobacterales',
'Gammaproteobacteria',
'Proteobacteria',
'Bacteria'],
't_journal': ['Acta Crystallogr. D Biol. Crystallogr.'],
't_mesh_terms': ['Crystallography, X-Ray,Humans,Hydrophobic and Hydrophilic '
'Interactions,Models, Molecular,NADP,Protein Structure, '
'Tertiary,Pyrimidines,Stereoisomerism,Tetrahydrofolate '
'Dehydrogenase'],
't_molecule_info': ['protein structure',
'homo',
'monomer',
'Dihydrofolate reductase',
'Dihydrofolate reductase',
'protein structure',
'homo',
'monomer',
'DHFR',
'P00374',
'P00374-2',
'Protein has possible alternate isoforms',
'DHFR',
'Protein has possible natural variant ',
'Nucleotide binding - NADP',
'DYR_HUMAN',
'DYR_HUMAN'],
't_molecule_sequence': 'VGSLNCIVAVSQNMGIGKNGDLPWPPLRNEFRYFQRMTTTSSVEGKQNLVIMGKKTWFSIPEKNRPLKGRINLVLSRELKEPPQGAHFLSRSLDDALKLTEQPELANKVDMVWIVGGSSVYKEAMNHPGHLKLFVTRIMQDFESDTFFPEIDLEKYKLLPEYPGVLSDVQEEKGIKYKFEVYEKND',
't_organism_name': ['Homo sapiens',
'Man',
'Homo Sapiens (Human)',
'Human',
'Homo Sapiens',
'Homo sapiens',
'Homo',
'Homininae',
'Hominidae',
'Primates',
'Mammalia',
'Chordata',
'Metazoa',
'Eukaryota'],
'tax_id': [9606],
'tax_query': [9606],
'title': 'Preferential Selection of Isomer Binding from Chiral Mixtures: '
'Alternate Binding Modes Observed for the E- and Z-isomers of a '
'Series of 5-Substituted 2,4-Diaminofuro[2,3-d]pyrimidines as '
'Ternary Complexes with NADPH and Human Dihydrofolate Reductase',
'uniprot': ['P00374 : DYR_HUMAN'],
'uniprot_accession': ['P00374', 'P00374-2'],
'uniprot_accession_best': ['P00374'],
'uniprot_best': ['P00374 : DYR_HUMAN'],
'uniprot_coverage': [0.99],
'uniprot_features': ['Protein has possible alternate isoforms',
'DHFR',
'Protein has possible natural variant ',
'Nucleotide binding - NADP'],
'uniprot_id': ['DYR_HUMAN', 'DYR_HUMAN'],
'uniprot_id_best': ['DYR_HUMAN'],
'uniprot_non_canonical': ['P00374-2 : DYR_HUMAN'],
'unp_count': 1,
'unp_nf90_accession': ['B0YJ76',
'A0A2K6C3Y8',
'S5WD14',
'S5VM81',
'P00374',
'A0A024RAQ3'],
'unp_nf90_id': ['B0YJ76_HUMAN',
'A0A2K6C3Y8_MACNE',
'S5WD14_SHISS',
'S5VM81_ECO57',
'DYR_HUMAN',
'A0A024RAQ3_HUMAN'],
'unp_nf90_organism': ['Homo sapiens (Human)',
'Macaca nemestrina (Pig-tailed macaque)',
'Shigella sonnei (strain Ss046)',
'Escherichia coli O157:H7',
'Homo sapiens (Human)',
'Homo sapiens (Human)'],
'unp_nf90_protein_name': ['Dihydrofolate reductase',
'DHFR domain-containing protein',
'Dihydrofolate reductase',
'Dihydrofolate reductase',
'Dihydrofolate reductase',
'Dihydrofolate reductase'],
'unp_nf90_tax_id': ['9606', '9545', '300269', '83334', '9606', '9606']}
A full list of available search terms is available using the following command.
[170]:
print('There are {} available search terms'.format(len(second_results[0].keys())))
print(second_results[0].keys())
keys_without_q = [q for q in second_results[0].keys() if not (q.startswith('q_') or (q.startswith('t_')))]
print(len(keys_without_q))
There are 448 available search terms
dict_keys(['abstracttext_unassigned', 't_abstracttext_unassigned', 'q_abstracttext_unassigned', 'all_assembly_composition', 't_molecule_info', 'q_all_assembly_composition', 'all_assembly_form', 'q_all_assembly_form', 'all_assembly_id', 'q_all_assembly_id', 'all_assembly_mol_wt', 'q_all_assembly_mol_wt', 'all_assembly_type', 'q_all_assembly_type', 'all_authors', 'q_all_authors', 'all_molecule_names', 'q_all_molecule_names', 'all_num_interacting_entity_id', 'q_all_num_interacting_entity_id', 'assembly_composition', 'q_assembly_composition', 'assembly_form', 'q_assembly_form', 'assembly_id', 'q_assembly_id', 'assembly_mol_wt', 'q_assembly_mol_wt', 'assembly_num_component', 'q_assembly_num_component', 'assembly_type', 'q_assembly_type', 'beam_source_name', 'q_beam_source_name', 'biological_cell_component', 'all_go_terms', 'q_all_go_terms', 't_all_go_terms', 'q_biological_cell_component', 'biological_function', 'q_biological_function', 'biological_process', 'q_biological_process', 'bound_compound_id', 'q_bound_compound_id', 'bound_compound_name', 'q_bound_compound_name', 'bound_compound_synonym', 'q_bound_compound_synonym', 'bound_compound_systematic_name', 'q_bound_compound_systematic_name', 'bound_compound_weight', 'q_bound_compound_weight', 'cath_architecture', 'all_structure_family', 'q_all_structure_family', 't_all_structure_family', 'q_cath_architecture', 'cath_class', 'q_cath_class', 'cath_code', 'q_cath_code', 'cath_homologous_superfamily', 'q_cath_homologous_superfamily', 'cath_topology', 'q_cath_topology', 'cell_a', 'q_cell_a', 'cell_alpha', 'q_cell_alpha', 'cell_b', 'q_cell_b', 'cell_beta', 'q_cell_beta', 'cell_c', 'q_cell_c', 'cell_gamma', 'q_cell_gamma', 'chain_id', 'q_chain_id', 'citation_authors', 't_citation_authors', 'q_citation_authors', 'citation_doi', 'q_citation_doi', 'citation_title', 't_citation_title', 'q_citation_title', 'citation_year', 'q_citation_year', 'cofactor_class', 't_all_compound_names', 'q_cofactor_class', 'cofactor_id', 'q_cofactor_id', 'compound_id', 'q_compound_id', 'compound_name', 'all_compound_names', 'q_all_compound_names', 'q_compound_name', 'compound_synonym', 'q_compound_synonym', 'compound_systematic_name', 'q_compound_systematic_name', 'compound_weight', 'q_compound_weight', 'crystallisation_cond', 'q_crystallisation_cond', 'crystallisation_method', 'q_crystallisation_method', 'crystallisation_ph', 'q_crystallisation_ph', 'crystallisation_reservoir', 't_entry_info', 'q_crystallisation_reservoir', 'crystallisation_temperature', 'q_crystallisation_temperature', 'data_collection_date', 'q_data_collection_date', 'data_collection_year', 'q_data_collection_year', 'data_quality', 'q_data_quality', 'data_reduction_software', 'q_data_reduction_software', 'data_scaling_software', 'q_data_scaling_software', 'deposition_date', 'q_deposition_date', 'deposition_site', 'q_deposition_site', 'deposition_year', 'q_deposition_year', 'detector', 'q_detector', 'detector_type', 'q_detector_type', 'diffraction_protocol', 'q_diffraction_protocol', 'diffraction_source_type', 'q_diffraction_source_type', 'diffraction_wavelengths', 'q_diffraction_wavelengths', 'ec_hierarchy_name', 'all_enzyme_names', 'q_all_enzyme_names', 't_all_enzyme_names', 'q_ec_hierarchy_name', 'ec_number', 'q_ec_number', 'entity_id', 'q_entity_id', 'entity_weight', 'q_entity_weight', 'entry_author_list', 'q_entry_author_list', 'entry_authors', 't_entry_authors', 'q_entry_authors', 'entry_entity', 'entry_lig_entity', 'q_entry_lig_entity', 'entry_organism_scientific_name', 'entry_uniprot', 'q_entry_uniprot', 'entry_uniprot_accession', 'q_entry_uniprot_accession', 'entry_uniprot_id', 'q_entry_uniprot_id', 'enzyme_name', 'q_enzyme_name', 'enzyme_num_name', 'q_enzyme_num_name', 'enzyme_systematic_name', 'q_enzyme_systematic_name', 'experiment_data_available', 'q_experiment_data_available', 'experimental_method', 'q_experimental_method', 'expression_host_sci_name', 'expression_organism_name', 'q_expression_organism_name', 't_expression_organism_name', 'q_expression_host_sci_name', 'expression_host_synonyms', 'q_expression_host_synonyms', 'expression_host_tax_id', 'q_expression_host_tax_id', 'gene_name', 'q_gene_name', 'genus', 'q_genus', 'go_id', 'q_go_id', 'go_mapping', 'q_go_mapping', 'has_bound_molecule', 'q_has_bound_molecule', 'has_carb_polymer', 'q_has_carb_polymer', 'has_cofactor', 'q_has_cofactor', 'has_modified_residues', 'q_has_modified_residues', 'has_reactant', 'q_has_reactant', 'homologus_pdb_entity_id', 'q_homologus_pdb_entity_id', 'int_lig_entity', 'q_int_lig_entity', 'interacting_ligands', 'q_interacting_ligands', 'interpro', 'all_sequence_family', 'q_all_sequence_family', 't_all_sequence_family', 'q_interpro', 'interpro_accession', 'q_interpro_accession', 'interpro_dom', 'q_interpro_dom', 'interpro_dom_acc', 'q_interpro_dom_acc', 'interpro_dom_name', 'q_interpro_dom_name', 'interpro_name', 'q_interpro_name', 'interpro_supfam', 'q_interpro_supfam', 'interpro_supfam_acc', 'q_interpro_supfam_acc', 'interpro_supfam_name', 'q_interpro_supfam_name', 'inv_overall_quality', 'q_inv_overall_quality', 'journal', 't_journal', 'q_journal', 'journal_page', 'q_journal_page', 'journal_volume', 'q_journal_volume', 'matthews_coefficient', 'q_matthews_coefficient', 'max_observed_residues', 'q_max_observed_residues', 'mesh_terms', 't_mesh_terms', 'q_mesh_terms', 'model_quality', 'q_model_quality', 'modified_residue_flag', 'q_modified_residue_flag', 'molecule_name', 'q_molecule_name', 'molecule_sequence', 'q_molecule_sequence', 't_molecule_sequence', 'molecule_synonym', 'q_molecule_synonym', 'molecule_type', 'q_molecule_type', 'mutation', 'q_mutation', 'nigli_cell_a', 'q_nigli_cell_a', 'nigli_cell_alpha', 'q_nigli_cell_alpha', 'nigli_cell_b', 'q_nigli_cell_b', 'nigli_cell_beta', 'q_nigli_cell_beta', 'nigli_cell_c', 'q_nigli_cell_c', 'nigli_cell_gamma', 'q_nigli_cell_gamma', 'nigli_cell_symmetry', 'q_nigli_cell_symmetry', 'num_interacting_entity_id', 'q_num_interacting_entity_id', 'num_r_free_reflections', 'q_num_r_free_reflections', 'number_of_bound_entities', 'q_number_of_bound_entities', 'number_of_bound_molecules', 'q_number_of_bound_molecules', 'number_of_copies', 'q_number_of_copies', 'number_of_models', 'q_number_of_models', 'number_of_polymer_entities', 'q_number_of_polymer_entities', 'number_of_polymer_residues', 'q_number_of_polymer_residues', 'number_of_polymers', 'q_number_of_polymers', 'number_of_protein_chains', 'q_number_of_protein_chains', 'organism_scientific_name', 'organism_name', 'q_organism_name', 't_organism_name', 'q_organism_scientific_name', 'organism_synonyms', 'q_organism_synonyms', 'overall_quality', 'q_overall_quality', 'pdb_accession', 'q_pdb_accession', 'pdb_format_compatible', 'q_pdb_format_compatible', 'pdb_id', 'q_pdb_id', 'percent_solvent', 'q_percent_solvent', 'pfam', 'q_pfam', 'pfam_accession', 'q_pfam_accession', 'pfam_clan', 'q_pfam_clan', 'pfam_clan_name', 'q_pfam_clan_name', 'pfam_description', 'q_pfam_description', 'pfam_name', 'q_pfam_name', 'pivot_resolution', 'q_pivot_resolution', 'polymer_length', 'q_polymer_length', 'prefered_assembly_id', 'q_prefered_assembly_id', 'primary_wavelength', 'q_primary_wavelength', 'processing_site', 'q_processing_site', 'pubmed_author_list', 'q_pubmed_author_list', 'pubmed_authors', 'q_pubmed_authors', 'pubmed_id', 'q_pubmed_id', 'r_factor', 'q_r_factor', 'r_free', 'q_r_free', 'r_work', 'q_r_work', 'rank', 'q_rank', 'reactant_id', 'q_reactant_id', 'refinement_software', 'q_refinement_software', 'release_date', 'q_release_date', 'release_year', 'q_release_year', 'resolution', 'q_resolution', 'revision_date', 'q_revision_date', 'revision_year', 'q_revision_year', 'sample_preparation_method', 'q_sample_preparation_method', 'seq_100_cluster_number', 'q_seq_100_cluster_number', 'seq_100_cluster_rank', 'q_seq_100_cluster_rank', 'seq_30_cluster_number', 'q_seq_30_cluster_number', 'seq_30_cluster_rank', 'q_seq_30_cluster_rank', 'seq_40_cluster_number', 'q_seq_40_cluster_number', 'seq_40_cluster_rank', 'q_seq_40_cluster_rank', 'seq_50_cluster_number', 'q_seq_50_cluster_number', 'seq_50_cluster_rank', 'q_seq_50_cluster_rank', 'seq_70_cluster_number', 'q_seq_70_cluster_number', 'seq_70_cluster_rank', 'q_seq_70_cluster_rank', 'seq_90_cluster_number', 'q_seq_90_cluster_number', 'seq_90_cluster_rank', 'q_seq_90_cluster_rank', 'seq_95_cluster_number', 'q_seq_95_cluster_number', 'seq_95_cluster_rank', 'q_seq_95_cluster_rank', 'serial_crystal_experiment', 'q_serial_crystal_experiment', 'spacegroup', 'q_spacegroup', 'status', 'q_status', 'struct_asym_id', 'q_struct_asym_id', 'structure_determination_method', 'q_structure_determination_method', 'structure_solution_software', 'q_structure_solution_software', 'superkingdom', 'q_superkingdom', 'tax_id', 'q_tax_id', 'tax_query', 'q_tax_query', 'title', 't_entry_title', 'q_title', 'uniprot', 'q_uniprot', 'uniprot_accession', 'q_uniprot_accession', 'uniprot_accession_best', 'q_uniprot_accession_best', 'uniprot_best', 'q_uniprot_best', 'uniprot_coverage', 'q_uniprot_coverage', 'uniprot_features', 'q_uniprot_features', 'uniprot_id', 'q_uniprot_id', 'uniprot_id_best', 'q_uniprot_id_best', 'uniprot_non_canonical', 'q_uniprot_non_canonical', 'unp_nf90_accession', 'q_unp_nf90_accession', 'unp_nf90_id', 'q_unp_nf90_id', 'unp_nf90_organism', 'q_unp_nf90_organism', 'unp_nf90_protein_name', 'q_unp_nf90_protein_name', 'unp_nf90_tax_id', 'q_unp_nf90_tax_id', '_version_', 'q_unp_count', 'unp_count'])
217
As you can see we get lots of data back about the individual molecule we have searched for and the PDB entries in which it is contained.
We can get the PDB ID and experimental method for this first row as follows. Note that experimental method is a list
[171]:
print(second_results[0].get('pdb_id'))
print(second_results[0].get('experimental_method'))
3nxv
['X-ray diffraction']
We can restrict the results to only the information we want using a filter so its easier to see the information we want.
[172]:
print('3rd search')
search_terms = Q(molecule_name="Dihydrofolate reductase",organism_name="Human")
filter_terms = ['pdb_id', 'experimental_method']
third_results = run_search(search_terms, filter_terms)
pprint(third_results)
3rd search
Number of results: 79
[{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3nxv'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5hsu'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5ht4'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5ht5'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5hqz'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5hvb'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5hui'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5hpb'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4m6k'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4qjc'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4keb'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4m6j'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4g95'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4qhv'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4kfj'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4m6l'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4kak'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4kbn'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4kd7'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5hsr'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1kms'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1kmv'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1drf'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1pd9'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1dhf'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1hfr'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1hfq'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '4ddr'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5hve'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '5hqy'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1dls'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1dlr'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1hfp'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1ohk'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1pd8'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1mvs'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1ohj'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1s3w'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3f8y'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3ghc'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '2w3a'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '2w3m'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3s7a'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3f8z'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3s3v'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3nzd'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3nxo'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3n0h'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1u71'},
{'experimental_method': ['Solution NMR'], 'pdb_id': '1yho'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1mvt'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1boz'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1s3v'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1pdb'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1u72'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '1s3u'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '2dhf'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '2c2t'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '2c2s'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '2w3b'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3ghw'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3eig'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3ghv'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3gi2'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3nxy'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3oaf'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3fs6'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3nxt'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3ntz'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3nxr'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3nxx'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3l3r'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3gyf'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3f91'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '3nu0'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '6de4'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '6dav'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '6a7e'},
{'experimental_method': ['X-ray diffraction'], 'pdb_id': '6a7c'}]
6) Analysing and plotting the results¶
We are going to use a Python package called Pandas to help us sort and visualise the results
First we have to do a bit of housekeeping, some of the results are lists (a PDB entry can have more than one experimental method or organism for example) so we need to change them into strings so we can use them in a graph
[173]:
def change_lists_to_strings(results):
"""
input - list of results from search
output - list of results with lists changed into strings
"""
for row in results:
for data in row:
if type(row[data]) == list:
# if there are any numbers in the list change them into strings
row[data] = [str(a) for a in row[data]]
# unique and sort the list and then change the list into a string
row[data] = ','.join(sorted(list(set(row[data]))))
return results
[174]:
results = change_lists_to_strings(third_results)
pprint(results)
[{'experimental_method': 'X-ray diffraction', 'pdb_id': '3nxv'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5hsu'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5ht4'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5ht5'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5hqz'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5hvb'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5hui'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5hpb'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4m6k'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4qjc'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4keb'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4m6j'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4g95'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4qhv'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4kfj'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4m6l'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4kak'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4kbn'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4kd7'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5hsr'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1kms'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1kmv'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1drf'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1pd9'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1dhf'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1hfr'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1hfq'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '4ddr'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5hve'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '5hqy'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1dls'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1dlr'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1hfp'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1ohk'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1pd8'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1mvs'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1ohj'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1s3w'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3f8y'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3ghc'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '2w3a'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '2w3m'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3s7a'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3f8z'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3s3v'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3nzd'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3nxo'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3n0h'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1u71'},
{'experimental_method': 'Solution NMR', 'pdb_id': '1yho'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1mvt'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1boz'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1s3v'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1pdb'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1u72'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '1s3u'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '2dhf'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '2c2t'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '2c2s'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '2w3b'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3ghw'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3eig'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3ghv'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3gi2'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3nxy'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3oaf'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3fs6'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3nxt'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3ntz'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3nxr'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3nxx'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3l3r'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3gyf'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3f91'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '3nu0'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '6de4'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '6dav'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '6a7e'},
{'experimental_method': 'X-ray diffraction', 'pdb_id': '6a7c'}]
Notice that the only thing that changed is [‘X-ray diffraction’] is now ‘X-ray diffraction’
If we wanted to know the experimental methods used to determine structures of Human Dihydrofolate reductase we could loop through the results and count how many entries use each experimental method.
We can use a Python package called Pandas to do this for us. It changes the results into a mini database - called a DataFrame.
[175]:
def pandas_dataset(list_of_results):
results = change_lists_to_strings(list_of_results) # we have added our function to change lists to strings
df = pd.DataFrame(results)
return df
df = pandas_dataset(list_of_results=results)
print(df)
experimental_method pdb_id
0 X-ray diffraction 3nxv
1 X-ray diffraction 5hsu
2 X-ray diffraction 5ht4
3 X-ray diffraction 5ht5
4 X-ray diffraction 5hqz
.. ... ...
74 X-ray diffraction 3nu0
75 X-ray diffraction 6de4
76 X-ray diffraction 6dav
77 X-ray diffraction 6a7e
78 X-ray diffraction 6a7c
[79 rows x 2 columns]
We can use the this to count how many PDB codes there are for each experimental method This groups PDB IDs by experimental method and then counts the number of unique PDB IDs per method.
[176]:
ds = df.groupby('experimental_method')['pdb_id'].nunique()
print(ds)
experimental_method
Solution NMR 1
X-ray diffraction 78
Name: pdb_id, dtype: int64
We can find which experimental method has the greatest (max) or lowest (min) number of entries.
[177]:
dt = ds.max()
print(dt)
dt = ds.min()
print(dt)
78
1
We can sort the results so its in decending order and then the first value is the experimental method with the highest number of results
[178]:
ds.sort_values(ascending=False).index[0]
[178]:
'X-ray diffraction'
Or sort ascending so the experimental method with the lowest number of results is given
[179]:
ds.sort_values(ascending=True).index[0]
[179]:
'Solution NMR'
Or we can then very easily plot these results as a bar chart
[180]:
ds.iplot(kind='bar')
We will make this into two functions so we can use them again
[181]:
def pandas_count(list_of_results, column_to_group_by):
df = pandas_dataset(list_of_results)
ds = df.groupby(column_to_group_by)['pdb_id'].nunique()
return ds
def pandas_min_max(list_of_results, column_to_group_by, get_min=True):
df = pandas_dataset(list_of_results)
if get_min:
ds = df.groupby(column_to_group_by)['pdb_id'].min()
else:
ds = df.groupby(column_to_group_by)['pdb_id'].max()
return ds
def pandas_plot(list_of_results, column_to_group_by, graph_type='bar'):
ds = pandas_count(list_of_results=list_of_results, column_to_group_by=column_to_group_by)
ds.iplot(kind=graph_type)
One for counting the results
[182]:
pandas_count(list_of_results=results, column_to_group_by='experimental_method')
[182]:
experimental_method
Solution NMR 1
X-ray diffraction 78
Name: pdb_id, dtype: int64
One for getting min or max
[183]:
print('updated search')
search_terms = Q(molecule_name="Dihydrofolate reductase",organism_name="Human")
filter_terms = ['pdb_id', 'resolution']
new_results = run_search(search_terms, filter_terms)
pandas_min_max(list_of_results=new_results, column_to_group_by='resolution')
updated search
Number of results: 79
[183]:
resolution
1.050 1kmv
1.090 1kms
1.201 4m6j
1.210 5hsr
1.230 3fs6
1.240 3ghw
1.270 2w3b
1.300 3ghc
1.350 3ntz
1.396 4m6k
1.400 2c2s
1.450 3f8y
1.460 5hqy
1.500 2c2t
1.530 3gi2
1.550 6dav
1.600 2w3m
1.610 4qhv
1.620 4qjc
1.650 5hpb
1.700 3eig
1.760 4kfj
1.800 1mvt
1.840 4kbn
1.850 6a7e
1.900 1mvs
1.920 3n0h
2.000 1drf
2.010 3f8z
2.050 4ddr
2.060 6a7c
2.100 1boz
2.200 1pd9
2.300 1dhf
2.411 6de4
2.500 1ohj
2.715 4kd7
Name: pdb_id, dtype: object
and one for plotting the results
[184]:
pandas_plot(list_of_results=results, column_to_group_by='experimental_method')
Remember this only searched through the first 10 results. To increase the number of entries we have to run the search again, this time setting number_of_rows to a number in the function run_search.
[185]:
search_terms = Q(molecule_name="Dihydrofolate reductase",organism_name="Human")
results = run_search(search_terms, number_of_rows=10000)
Number of results: 79
Then we can count the results using our pandas function above
[186]:
pandas_count(list_of_results=results, column_to_group_by='experimental_method')
[186]:
experimental_method
Solution NMR 1
X-ray diffraction 78
Name: pdb_id, dtype: int64
Changing the result so it groups by release year of the PDB entries.
[187]:
pandas_count(list_of_results=results, column_to_group_by='release_year')
[187]:
release_year
1990 2
1992 1
1995 2
1998 6
2002 2
2003 5
2004 3
2005 3
2007 2
2009 13
2010 7
2011 7
2012 1
2013 8
2014 1
2015 2
2017 10
2018 2
2019 2
Name: pdb_id, dtype: int64
And then plot the number of entries released per year
[188]:
pandas_plot(list_of_results=results, column_to_group_by='release_year')
We can make this into a line graph
[189]:
pandas_plot(list_of_results=results, column_to_group_by='release_year', graph_type='line')
Try changing the term you want to search for and see if you get interesting results.
7) searching for two terms at once¶
It would be interesting to see how many PDB entries were solved by each experimental method per year.
we can use the tag “release_year” to get the year of release of each entry
We have to define a new function to group entries by two terms.
When we do the search we have to filter the results by the terms we want to plot otherwise it takes too long to run.
[190]:
search_terms = Q(all_enzyme_names="Lysozyme")
filter_results = ['beam_source_name','release_year', 'pdb_id']
results = run_search(search_terms, filter_results, number_of_rows=10000)
Number of results: 1975
This will take a while as it will return lots of results. We can then define a function to group the results by two terms.
[191]:
def pandas_plot_multi_groupby(results, first_column_to_group_by, second_column_to_group_by, y_axis='pdb_id', graph_type='line'):
df = pandas_dataset(results)
new_df = df.groupby([first_column_to_group_by, second_column_to_group_by])
ds = new_df.count().unstack().reset_index(first_column_to_group_by)
ds.plot(x=first_column_to_group_by, y=y_axis, kind=graph_type).legend(bbox_to_anchor=(1,1))
def pandas_plot_multi_groupby_min(results, first_column_to_group_by, second_column_to_group_by, graph_type='line', use_min=False, use_max=False):
df = pandas_dataset(results)
new_df = df.groupby([first_column_to_group_by])[second_column_to_group_by]
ds = None
if use_min:
ds = new_df.min()
elif use_max:
ds = new_df.max()
else:
print('specify either use_min or use_max')
return None
ds.plot(x=first_column_to_group_by, y=second_column_to_group_by, kind=graph_type)
def pandas_box_plot(results, first_column_to_group_by, second_column_to_group_by):
df = pandas_dataset(results)
df.boxplot(column=second_column_to_group_by,by=first_column_to_group_by)
[192]:
pandas_plot_multi_groupby(results, 'release_year', 'beam_source_name')
This shows us that rotating anodes were used as the major source of radiation until around 2004 when Synchrotron’s overtook as the major source of radiation.
Try editing the queries to plot interesting trends within the PDB
Questions to answer¶
Electron Microscopy is going through a revolution. Is this leading to a growth in Electron Microscopy PDB entries?
[193]:
search_terms = Q(experimental_method='X-ray diffraction',release_year='2018')
filter_results = ['structure_determination_method','release_year', 'pdb_id']
results = run_search(search_terms, filter_results, number_of_rows=1000000)
pandas_plot(list_of_results=results, column_to_group_by='structure_determination_method')
Number of results: 15794
New refinement programs have got better and there are more methods to validate the quality of strucures in the PDB. Have structures got better over time? We can use “overall_quality” to judge this This could be plotted as a groupby or a box plot.
[194]:
search_terms = Q(experimental_method='Electron Microscopy')
filter_results = ['experimental_method','release_year', 'pdb_id']
results = run_search(search_terms, filter_results, number_of_rows=1000000)
pandas_plot_multi_groupby(results, 'release_year', 'experimental_method')
Number of results: 61831
Electron Microscopy resolution has been said to be improving. Is this true? hint - the search term and filter can be different. pandas_plot_multi_groupby_min with use_min would be useful to plot this or maybe a box plot?
[195]:
search_terms = Q(experimental_method='Electron Microscopy')
filter_results = ['overall_quality','release_year', 'pdb_id']
results = run_search(search_terms, filter_results, number_of_rows=1000000)
pandas_plot_multi_groupby_min(results, 'release_year', 'overall_quality', use_min=True)
pandas_box_plot(results, 'release_year', 'overall_quality')
Number of results: 61831
It has been said that all the simple structures have been done and that only complicated structures are left. One metric for “complicated” could be size of the assembly - assembly_mol_wt would be useful here.
[196]:
search_terms = Q(experimental_method='Electron Microscopy')
filter_results = ['assembly_mol_wt','release_year', 'pdb_id']
results = run_search(search_terms, filter_results, number_of_rows=1000000)
pandas_box_plot(results, 'release_year', 'assembly_mol_wt')
Number of results: 61831