The Protein Data Bank (PDB) is the primary global archive for experimentally determined three-dimensional structures of biological macromolecules. The RCSB PDB exposes these data through programmatic interfaces that support search, metadata retrieval, coordinate download, and access to assembly- and entity-level annotations. For structural bioinformatics, these APIs make it possible to move from a biological question to a reproducible computational workflow without manual browsing of the PDB website.
rPDBapi provides an R interface to these services. It
combines search helpers, operator-based query construction, metadata
retrieval, GraphQL-based data fetching, and structure download utilities
in a form that fits naturally into R-based data analysis pipelines.
This vignette is written for users who want to retrieve and analyze PDB data directly in R. The examples focus on protein kinase structures because kinases are biologically important, structurally diverse, and common targets in drug-discovery workflows.
if (!live_examples) {
cat(
"This vignette contains live API examples. ",
"To execute them during rendering, set `RPDBAPI_RUN_LIVE=true` before knitting. ",
"When the flag is not set, the code is still shown but network-dependent chunks are not evaluated.\n",
sep = ""
)
}
This vignette contains live API examples. To execute them during
rendering, set RPDBAPI_RUN_LIVE=true before knitting. When
the flag is not set, the code is still shown but network-dependent
chunks are not evaluated.
The note above makes the vignette safe to build in offline or CRAN-like environments while keeping all examples directly executable in a live session. This is especially important for API-driven packages, where remote services are part of the normal workflow but should not make documentation brittle.
install.packages("rPDBapi")
# Development version
remotes::install_github("selcukorkmaz/rPDBapi")
The package can be installed from CRAN or from the development repository. The development version is useful when you want the newest API mappings, tests, and documentation updates. In a package vignette, installation code is normally shown but not run, because users execute it once interactively rather than during every render.
library(rPDBapi)
library(dplyr)
This chunk loads the package and dplyr, which we will
use for simple tabulation and ranking. Most structural bioinformatics
workflows combine API access with data manipulation, so it is useful to
establish that pattern from the beginning.
Programmatic PDB access is valuable when you need to:
In practice, this means that a question such as “Which high-resolution protein kinase structures are available, what organisms do they come from, and what do their assemblies look like?” can be answered in one analysis script instead of through manual web browsing.
At a high level, the package supports seven related tasks:
The current package version also hardens return contracts and error
classes. That matters when rPDBapi is embedded in larger
pipelines, because downstream code can now make stronger assumptions
about object classes, identifier formats, and failure modes.
rPDBapi exports functions that fall into nine practical
groups:
query_search(),
perform_search()DefaultOperator(),
ExactMatchOperator(), InOperator(),
ContainsWordsOperator(),
ContainsPhraseOperator(),
ComparisonOperator(), RangeOperator(),
ExistsOperator(), SequenceOperator(),
SeqMotifOperator(), StructureOperator(),
ChemicalOperator()QueryNode(),
QueryGroup(), RequestOptions(),
infer_search_service(), ScoredResult()infer_id_type(),
parse_rcsb_id(), build_entry_id(),
build_assembly_id(), build_entity_id(),
build_instance_id()data_fetcher(),
fetch_data(), generate_json_query(),
get_info(), find_results(),
find_papers(), describe_chemical(),
get_fasta_from_rcsb_entry()list_rcsb_fields(),
search_rcsb_fields(), validate_properties(),
add_property()data_fetcher_batch(), cache_info(),
clear_rpdbapi_cache()get_pdb_file(),
get_pdb_api_url()as_rpdb_entry(),
as_rpdb_assembly(), as_rpdb_polymer_entity(),
as_rpdb_chemical_component(),
as_rpdb_structure(), summarize_entries(),
summarize_assemblies(),
extract_taxonomy_table(),
extract_ligand_table(),
extract_calpha_coordinates(),
join_structure_sequence()send_api_request(), handle_api_errors(),
parse_response(), search_graphql(),
return_data_as_dataframe()The main workflow only needs a subset of these functions, but the full package is designed as a layered interface. High-level helpers are convenient for routine work, while low-level helpers make it possible to debug requests, build custom workflows, or extend the package into larger pipelines.
Before starting with code, it helps to distinguish a few PDB concepts:
ENTRY: a PDB deposition such as 4HHBASSEMBLY: a biological assembly within an entry, such
as 4HHB-1POLYMER_ENTITY: a unique macromolecular entity within
an entry, such as a protein chain definitionATTRIBUTE: a searchable or retrievable field, such as
exptl.method or
rcsb_entry_info.molecular_weightrPDBapi mirrors these levels. Search functions return
identifiers at the appropriate level, and metadata functions use those
identifiers to fetch the corresponding records.
kinase_full_text <- DefaultOperator("protein kinase")
high_resolution <- RangeOperator(
attribute = "rcsb_entry_info.resolution_combined",
from_value = 0,
to_value = 2.5
)
xray_method <- ExactMatchOperator(
attribute = "exptl.method",
value = "X-RAY DIFFRACTION"
)
kinase_query <- QueryGroup(
queries = list(kinase_full_text, xray_method, high_resolution),
logical_operator = "AND"
)
kinase_query
#> $type
#> [1] "group"
#>
#> $logical_operator
#> [1] "and"
#>
#> $nodes
#> $nodes[[1]]
#> $nodes[[1]]$type
#> [1] "terminal"
#>
#> $nodes[[1]]$service
#> [1] "full_text"
#>
#> $nodes[[1]]$parameters
#> $value
#> [1] "protein kinase"
#>
#> attr(,"class")
#> [1] "DefaultOperator" "list"
#>
#>
#> $nodes[[2]]
#> $nodes[[2]]$type
#> [1] "terminal"
#>
#> $nodes[[2]]$service
#> [1] "text"
#>
#> $nodes[[2]]$parameters
#> $attribute
#> [1] "exptl.method"
#>
#> $value
#> [1] "X-RAY DIFFRACTION"
#>
#> $operator
#> [1] "exact_match"
#>
#> attr(,"class")
#> [1] "ExactMatchOperator" "list"
#>
#>
#> $nodes[[3]]
#> $nodes[[3]]$type
#> [1] "terminal"
#>
#> $nodes[[3]]$service
#> [1] "text"
#>
#> $nodes[[3]]$parameters
#> $operator
#> [1] "range"
#>
#> $attribute
#> [1] "rcsb_entry_info.resolution_combined"
#>
#> $negation
#> [1] FALSE
#>
#> $value
#> $value$from
#> [1] 0
#>
#> $value$to
#> [1] 2.5
#>
#>
#> attr(,"class")
#> [1] "RangeOperator" "list"
This code builds a structured query object without contacting the API. The query says: search for records related to “protein kinase”, require X-ray diffraction as the experimental method, and restrict the results to structures with a reported resolution of 2.5 angstroms or better. This is a useful pattern because it separates biological intent from the mechanics of the HTTP request.
search_controls <- RequestOptions(
result_start_index = 0,
num_results = 10,
sort_by = "score",
desc = TRUE
)
search_controls
#> $paginate
#> $paginate$start
#> [1] 0
#>
#> $paginate$rows
#> [1] 10
#>
#>
#> $sort
#> $sort[[1]]
#> $sort[[1]]$sort_by
#> [1] "score"
#>
#> $sort[[1]]$direction
#> [1] "desc"
RequestOptions() defines how many hits to return and how
to sort them. In other words, the query object describes what you want,
and the request options describe how you want it delivered. That
distinction matters when you are iterating over result pages or creating
reproducible subsets for downstream analysis.
example_ids <- c("4HHB", "4HHB-1", "4HHB_1", "4HHB.A", "ATP")
dplyr::tibble(
id = example_ids,
inferred_type = infer_id_type(example_ids)
)
#> # A tibble: 5 × 2
#> id inferred_type
#> <chr> <chr>
#> 1 4HHB ENTRY
#> 2 4HHB-1 ASSEMBLY
#> 3 4HHB_1 ENTITY
#> 4 4HHB.A INSTANCE
#> 5 ATP CHEMICAL_COMPONENT
parse_rcsb_id("4HHB-1")
#> $raw_id
#> [1] "4HHB-1"
#>
#> $normalized_id
#> [1] "4HHB-1"
#>
#> $id_type
#> [1] "ASSEMBLY"
#>
#> $entry_id
#> [1] "4HHB"
#>
#> $assembly_id
#> [1] "1"
#>
#> $entity_id
#> NULL
#>
#> $instance_id
#> NULL
#>
#> $separator
#> [1] "-"
build_entry_id(" 4HHB ")
#> [1] "4HHB"
build_assembly_id("4HHB", 1)
#> [1] "4HHB-1"
build_entity_id("4HHB", 1)
#> [1] "4HHB_1"
build_instance_id("4HHB", "A")
#> [1] "4HHB.A"
These helpers make identifier handling explicit. They are useful when a workflow moves across entry-, assembly-, entity-, and chain-level retrieval, because the required identifier syntax changes with the biological level. In practice, the helpers reduce ad hoc string handling and make it easier to write validation checks before a request is sent.
kinase_hits <- query_search("protein kinase")
head(kinase_hits, 10)
class(kinase_hits)
attr(kinase_hits, "return_type")
query_search() is the fastest way to ask a general
question of the archive. Here it performs a full-text search and returns
entry identifiers. The returned object is not just a plain character
vector: it carries the class rPDBapi_query_ids, which makes
the contract explicit and helps downstream code reason about what kind
of object was returned.
kinase_entry_ids <- perform_search(
search_operator = kinase_query,
return_type = "ENTRY",
request_options = search_controls,
verbosity = FALSE
)
kinase_entry_ids
class(kinase_entry_ids)
perform_search() executes the operator-based query
assembled earlier. This is the function to use when you need precise
control over attributes, logical combinations, return types, or
pagination. In structural bioinformatics, this kind of targeted search
is often more useful than full-text search alone, because it lets you
combine biological meaning with experimental constraints.
entry_properties <- list(
rcsb_id = list(),
struct = c("title"),
struct_keywords = c("pdbx_keywords"),
exptl = c("method"),
rcsb_entry_info = c("molecular_weight", "resolution_combined"),
rcsb_accession_info = c("initial_release_date")
)
entry_properties
#> $rcsb_id
#> list()
#>
#> $struct
#> [1] "title"
#>
#> $struct_keywords
#> [1] "pdbx_keywords"
#>
#> $exptl
#> [1] "method"
#>
#> $rcsb_entry_info
#> [1] "molecular_weight" "resolution_combined"
#>
#> $rcsb_accession_info
#> [1] "initial_release_date"
This property list defines the fields we want from the GraphQL endpoint. It captures both structural metadata and biologically meaningful annotations: structure title, keywords, experimental method, molecular weight, resolution, and release date. By stating these fields explicitly, the workflow remains transparent and easy to reproduce.
head(list_rcsb_fields("ENTRY"), 10)
#> data_type field subfield
#> 1 ENTRY rcsb_id <NA>
#> 2 ENTRY struct title
#> 3 ENTRY struct_keywords pdbx_keywords
#> 4 ENTRY exptl method
#> 5 ENTRY cell length_a
#> 6 ENTRY cell length_b
#> 7 ENTRY cell length_c
#> 8 ENTRY cell volume
#> 9 ENTRY cell angle_beta
#> 10 ENTRY citation title
search_rcsb_fields("resolution", data_type = "ENTRY")
#> data_type field subfield
#> 13 ENTRY rcsb_entry_info resolution_combined
validate_properties(
properties = entry_properties,
data_type = "ENTRY",
strict = TRUE
)
validate_properties(
properties = list(
rcsb_entry_info = c("resolution_combined", "unknown_subfield")
),
data_type = "ENTRY",
strict = FALSE
)
#> $unknown_fields
#> character(0)
#>
#> $unknown_subfields
#> $unknown_subfields$rcsb_entry_info
#> [1] "unknown_subfield"
The schema-aware helpers are useful when building property lists
iteratively. list_rcsb_fields() exposes the package’s
built-in field registry, search_rcsb_fields() narrows it to
a topic of interest, and validate_properties() checks that
a property list matches the expected data-type-specific structure. In
strict mode, validation fails early; in non-strict mode, it returns
diagnostics that can be incorporated into an interactive workflow or a
package test.
old_opt <- options(rPDBapi.strict_property_validation = TRUE)
on.exit(options(old_opt), add = TRUE)
generate_json_query(
ids = c("4HHB"),
data_type = "ENTRY",
properties = list(rcsb_entry_info = c("resolution_combined"))
)
This option-gated strict mode is useful when you want a pipeline to reject unknown fields immediately. Because the option is opt-in, the package preserves backward compatibility for existing code while still supporting stricter validation in controlled workflows.
kinase_metadata <- data_fetcher(
id = kinase_entry_ids[1:5],
data_type = "ENTRY",
properties = entry_properties,
return_as_dataframe = TRUE
)
kinase_metadata
data_fetcher() is the main high-level retrieval function
for metadata. It takes identifiers, the data level of interest, and a
property list, then returns either a validated nested response or a
flattened data frame. For many analysis tasks, returning a data frame is
the most convenient choice because it fits directly into standard R
workflows for filtering, joining, and plotting.
kinase_json_query <- generate_json_query(
ids = kinase_entry_ids[1:3],
data_type = "ENTRY",
properties = entry_properties
)
cat(kinase_json_query)
This chunk exposes the GraphQL query string that rPDBapi
sends to the RCSB data API. Seeing the generated query is helpful when
you are debugging a field name, comparing package output with the
official schema, or teaching others how the package maps R objects to
API requests.
kinase_raw <- fetch_data(
json_query = kinase_json_query,
data_type = "ENTRY",
ids = kinase_entry_ids[1:3]
)
str(kinase_raw, max.level = 2)
fetch_data() returns a validated raw payload and tags it
with the class rPDBapi_fetch_response. This is useful when
you want to inspect nested JSON content before flattening it, preserve
hierarchy for custom parsing, or verify that a field is present before
building a larger workflow around it.
kinase_tidy <- return_data_as_dataframe(
response = kinase_raw,
data_type = "ENTRY",
ids = kinase_entry_ids[1:3]
)
kinase_tidy
return_data_as_dataframe() converts the nested response
into a rectangular R data structure. This transformation is central to
reproducible bioinformatics: once the results are tidy, they can be
analyzed with dplyr, joined to other annotations,
summarized statistically, or passed to visualization packages.
High-throughput structural workflows rarely stop at one or two identifiers. In screening, comparative analysis, or annotation projects, it is common to fetch dozens or hundreds of records with the same property specification.
cache_dir <- file.path(tempdir(), "rpdbapi-vignette-cache")
kinase_batch <- data_fetcher_batch(
id = kinase_entry_ids[1:5],
data_type = "ENTRY",
properties = entry_properties,
return_as_dataframe = TRUE,
batch_size = 2,
retry_attempts = 2,
retry_backoff = 0,
cache = TRUE,
cache_dir = cache_dir,
progress = FALSE,
verbosity = FALSE
)
kinase_batch
attr(kinase_batch, "provenance")
cache_info(cache_dir)
data_fetcher_batch() scales the single-request
data_fetcher() model to larger identifier sets. It splits
requests into batches, retries transient failures, optionally stores
batch results on disk, and attaches provenance to the returned object.
That provenance is important for reproducibility because it records the
retrieval mode, batch size, retry configuration, and cache usage.
clear_rpdbapi_cache(cache_dir)
cache_info(cache_dir)
This cache-management pattern is especially useful in iterative analysis. Repeated metadata retrieval becomes faster when the same requests are reused, while explicit cache inspection and cleanup keep the workflow transparent.
# Use data_fetcher() when:
# - the ID set is small
# - you want the simplest request path
# - retry, cache, and provenance are unnecessary
# Use data_fetcher_batch() when:
# - the ID set is large
# - requests may need retries
# - repeated retrieval should reuse cached results
# - you want an explicit provenance record
In practice, data_fetcher() is usually sufficient for
exploratory work. data_fetcher_batch() becomes more useful
as the workflow moves toward larger or repeated retrieval, where retry
behavior, caching, and provenance become part of the analysis design
rather than implementation detail.
provenance_tbl <- dplyr::tibble(
field = names(attr(kinase_batch, "provenance")),
value = vapply(
attr(kinase_batch, "provenance"),
function(x) {
if (is.list(x)) "<list>" else as.character(x)
},
character(1)
)
)
provenance_tbl
Interpreting provenance explicitly is useful when results are produced in a non-interactive workflow. The provenance record makes it clear how many batches were used, whether caching was enabled, and how the retrieval was configured, which makes the metadata table easier to audit later.
Biological assemblies are often the correct unit of interpretation for oligomeric proteins. A deposited asymmetric unit may not reflect the functional quaternary structure, so assembly-level retrieval is important when studying stoichiometry, symmetry, and interfaces.
kinase_assembly_ids <- perform_search(
search_operator = kinase_query,
return_type = "ASSEMBLY",
request_options = RequestOptions(result_start_index = 0, num_results = 5),
verbosity = FALSE
)
kinase_assembly_ids
This search requests assembly identifiers rather than entry identifiers. The returned IDs encode both the entry and the assembly number, making them appropriate inputs for assembly-level metadata retrieval. This is an important distinction in structural biology because entry-level and assembly-level questions are not interchangeable.
assembly_properties <- list(
rcsb_id = list(),
pdbx_struct_assembly = c("details", "method_details", "oligomeric_count"),
rcsb_struct_symmetry = c("kind", "symbol")
)
kinase_assemblies <- data_fetcher(
id = kinase_assembly_ids,
data_type = "ASSEMBLY",
properties = assembly_properties,
return_as_dataframe = TRUE
)
kinase_assemblies
This chunk retrieves assembly descriptors and symmetry annotations. In practice, these fields help answer questions about oligomeric state and biological interpretation, such as whether a kinase structure is monomeric, dimeric, or associated with a symmetric assembly.
assembly_object <- as_rpdb_assembly(
kinase_assemblies,
metadata = list(query = "protein kinase assemblies")
)
assembly_object
dplyr::as_tibble(assembly_object)
summarize_assemblies(assembly_object)
The assembly object wrapper is useful when you want to retain
lightweight metadata alongside a table while still working with
tibble-oriented tools. summarize_assemblies() then provides
a narrow helper for common assembly questions, such as typical
oligomeric count and the diversity of symmetry labels in the retrieved
result set.
One practical source of bugs in structural workflows is mixing
identifier levels. A valid entry ID is not automatically a valid
assembly or entity ID, and the corresponding data_type must
match the biological level of the request.
dplyr::tibble(
example_id = c("4HHB", "4HHB-1", "4HHB_1", "4HHB.A", "ATP"),
inferred_type = infer_id_type(c("4HHB", "4HHB-1", "4HHB_1", "4HHB.A", "ATP"))
)
#> # A tibble: 5 × 2
#> example_id inferred_type
#> <chr> <chr>
#> 1 4HHB ENTRY
#> 2 4HHB-1 ASSEMBLY
#> 3 4HHB_1 ENTITY
#> 4 4HHB.A INSTANCE
#> 5 ATP CHEMICAL_COMPONENT
parse_rcsb_id("4HHB.A")
#> $raw_id
#> [1] "4HHB.A"
#>
#> $normalized_id
#> [1] "4HHB.A"
#>
#> $id_type
#> [1] "INSTANCE"
#>
#> $entry_id
#> [1] "4HHB"
#>
#> $assembly_id
#> NULL
#>
#> $entity_id
#> NULL
#>
#> $instance_id
#> [1] "A"
#>
#> $separator
#> [1] "."
This is useful as a preflight step before retrieval. In larger workflows, a small identifier check often saves time because it catches level mismatches before the request reaches the API.
# Entry-level retrieval
data_fetcher(
id = build_entry_id("4HHB"),
data_type = "ENTRY",
properties = list(rcsb_id = list())
)
# Assembly-level retrieval
data_fetcher(
id = build_assembly_id("4HHB", 1),
data_type = "ASSEMBLY",
properties = list(rcsb_id = list())
)
# Polymer-entity retrieval
data_fetcher(
id = build_entity_id("4HHB", 1),
data_type = "POLYMER_ENTITY",
properties = list(rcsb_id = list())
)
Using the builder helpers makes the intended record level explicit in code. That is especially helpful when identifiers are generated programmatically from entry IDs and entity or assembly indices.
Many analyses need entity-level annotations rather than whole-entry metadata. For example, taxonomy belongs naturally to polymer entities because different entities within the same structure can come from different organisms or constructs.
kinase_polymer_ids <- perform_search(
search_operator = kinase_query,
return_type = "POLYMER_ENTITY",
request_options = RequestOptions(result_start_index = 0, num_results = 5),
verbosity = FALSE
)
kinase_polymer_ids
Here the same biological query is projected onto polymer entities. This is useful when you want annotations at the chain-definition level, such as source organism, sequence grouping, or entity-specific descriptors.
polymer_properties <- list(
rcsb_id = list(),
rcsb_entity_source_organism = c("ncbi_taxonomy_id", "ncbi_scientific_name"),
rcsb_cluster_membership = c("cluster_id", "identity")
)
kinase_polymer_metadata <- data_fetcher(
id = kinase_polymer_ids,
data_type = "POLYMER_ENTITY",
properties = polymer_properties,
return_as_dataframe = TRUE
)
kinase_polymer_metadata
This result provides organism-level context that is often essential in comparative structural biology. For example, you might use these fields to separate human kinase structures from bacterial homologs, or to identify closely related entities before selecting representatives for downstream modeling.
polymer_object <- as_rpdb_polymer_entity(
kinase_polymer_metadata,
metadata = list(query = "kinase polymer entities")
)
taxonomy_table <- extract_taxonomy_table(polymer_object)
taxonomy_table
taxonomy_table %>%
count(ncbi_scientific_name, sort = TRUE)
extract_taxonomy_table() is intentionally narrow: it
keeps only the fields needed to represent source-organism assignments
cleanly. This is useful when a larger polymer-entity table contains many
retrieval columns, but the immediate analysis question is taxonomic
composition or species-level redundancy.
selected_entry <- kinase_entry_ids[[1]]
selected_info <- get_info(selected_entry)
entry_summary <- dplyr::tibble(
rcsb_id = selected_entry,
title = purrr::pluck(selected_info, "struct", "title", .default = NA_character_),
keywords = purrr::pluck(selected_info, "struct_keywords", "pdbx_keywords", .default = NA_character_),
method = purrr::pluck(selected_info, "exptl", 1, "method", .default = NA_character_),
citation_title = purrr::pluck(selected_info, "rcsb_primary_citation", "title", .default = NA_character_),
resolution = paste(
purrr::pluck(selected_info, "rcsb_entry_info", "resolution_combined", .default = NA),
collapse = "; "
)
)
entry_summary
get_info() retrieves a full entry record as a nested
list. This is useful when you want a richer, less filtered
representation than a GraphQL property subset. In this example, we
extract structure title, keywords, experimental method, citation title,
and resolution to build a compact summary of one kinase entry. These
fields are exactly the kinds of annotations structural biologists
inspect when deciding whether a structure is suitable for biological
interpretation or downstream modeling.
kinase_papers <- find_papers("protein kinase", max_results = 3)
kinase_keywords <- find_results("protein kinase", field = "struct_keywords")
kinase_papers
utils::head(kinase_keywords, 3)
These helper functions show how rPDBapi can bridge
structures and biological interpretation. find_papers()
provides publication titles associated with matching entries, while
find_results() can retrieve selected metadata fields across
search results. Together, they support a workflow in which structures,
annotations, and literature references are linked inside the same R
session.
kinase_structure <- get_pdb_file(
pdb_id = selected_entry,
filetype = "cif",
verbosity = FALSE
)
coordinate_matrix <- matrix(kinase_structure$xyz, ncol = 3, byrow = TRUE)
coordinate_df <- data.frame(
x = coordinate_matrix[, 1],
y = coordinate_matrix[, 2],
z = coordinate_matrix[, 3]
)
calpha_atoms <- cbind(
kinase_structure$atom[kinase_structure$calpha, c("chain", "resno", "resid")],
coordinate_df[kinase_structure$calpha, , drop = FALSE]
)
head(calpha_atoms, 10)
get_pdb_file() downloads and parses the structure file
into an object that contains atomic records and coordinates. This is the
transition point from metadata analysis to coordinate analysis. The
example extracts C-alpha atoms, which are commonly used in structural
alignment, geometry summaries, coarse distance analyses, and quick
visual inspection of protein backbones.
calpha_atoms <- extract_calpha_coordinates(kinase_structure)
head(calpha_atoms, 10)
extract_calpha_coordinates() packages a common
structural-bioinformatics step into a reusable helper. The result is
immediately usable for plotting, distance-based summaries, or
chain-level coordinate analyses without manually reconstructing the
atom/coordinate join each time.
kinase_sequences <- get_fasta_from_rcsb_entry(selected_entry, verbosity = FALSE)
length(kinase_sequences)
utils::head(nchar(unlist(kinase_sequences)))
The FASTA workflow complements coordinate retrieval by exposing the underlying macromolecular sequences. Having both sequence and structure available in the same environment is useful for tasks such as domain boundary checks, sequence length summaries, and linking structural hits to sequence-based pipelines.
chain_sequence_summary <- join_structure_sequence(
kinase_structure,
kinase_sequences
)
chain_sequence_summary
This helper joins sequence-level and coordinate-level summaries at the chain level. In practice, it provides a quick diagnostic for whether the downloaded structure and the FASTA content align as expected, and it creates a compact table that can be extended with additional chain annotations.
The typed object wrappers are intentionally lightweight. They do not replace the underlying data; instead, they add a stable class layer for printing, conversion, and helper dispatch.
entry_demo <- as_rpdb_entry(
data.frame(
rcsb_id = c("4HHB", "1CRN"),
method = c("X-RAY DIFFRACTION", "SOLUTION NMR"),
resolution_combined = c("1.74", NA),
stringsAsFactors = FALSE
),
metadata = list(example = "local object demo")
)
entry_demo
#> <rPDBapi_entry> with data class: data.frame
dplyr::as_tibble(entry_demo)
#> # A tibble: 2 × 3
#> rcsb_id method resolution_combined
#> <chr> <chr> <chr>
#> 1 4HHB X-RAY DIFFRACTION 1.74
#> 2 1CRN SOLUTION NMR <NA>
summarize_entries(entry_demo)
#> # A tibble: 1 × 4
#> n_entries n_methods best_resolution median_molecular_weight
#> <int> <int> <dbl> <dbl>
#> 1 2 2 1.74 NA
entry_demo$metadata
#> $example
#> [1] "local object demo"
This pattern is useful when a workflow needs to preserve both data
and context. For example, metadata can record which query produced a
table or which processing choices were applied. The
as_tibble() methods then let the object drop back into a
standard tidyverse pipeline without extra conversion code.
structure_demo <- as_rpdb_structure(
list(
atom = data.frame(
chain = c("A", "A"),
resno = c(1L, 2L),
resid = c("GLY", "ALA"),
stringsAsFactors = FALSE
),
xyz = c(1, 2, 3, 4, 5, 6),
calpha = c(TRUE, FALSE)
),
metadata = list(source = "illustration")
)
structure_demo
#> <rPDBapi_structure> with data class: list
dplyr::as_tibble(structure_demo)
#> # A tibble: 6 × 3
#> atom xyz calpha
#> <chr> <chr> <chr>
#> 1 A;A;1;2;GLY;ALA 1 TRUE
#> 2 A;A;1;2;GLY;ALA 2 FALSE
#> 3 A;A;1;2;GLY;ALA 3 TRUE
#> 4 A;A;1;2;GLY;ALA 4 FALSE
#> 5 A;A;1;2;GLY;ALA 5 TRUE
#> 6 A;A;1;2;GLY;ALA 6 FALSE
The structure wrapper is especially useful when one analysis session contains multiple parsed structures and derived tables. The class layer makes those objects easier to distinguish and easier to handle consistently.
entry_object <- as_rpdb_entry(
kinase_metadata,
metadata = list(query = "protein kinase entry metadata")
)
summarize_entries(entry_object)
kinase_summary <- dplyr::as_tibble(entry_object) %>%
mutate(
molecular_weight = as.numeric(molecular_weight),
resolution_combined = as.numeric(resolution_combined),
initial_release_date = as.Date(initial_release_date)
) %>%
arrange(resolution_combined) %>%
select(
rcsb_id,
title,
pdbx_keywords,
method,
molecular_weight,
resolution_combined,
initial_release_date
)
kinase_summary
kinase_summary %>%
summarise(
n_structures = n(),
median_molecular_weight = median(molecular_weight, na.rm = TRUE),
best_resolution = min(resolution_combined, na.rm = TRUE)
)
Once the metadata are in a data frame, standard R analysis becomes
immediate. This chunk ranks the retrieved kinase structures by
resolution and computes a few simple summaries. Although the analysis is
straightforward, it illustrates the main advantage of using
rPDBapi: PDB metadata become ordinary tabular R data that
can be manipulated with the same tools used elsewhere in
bioinformatics.
kinase_polymer_metadata %>%
count(ncbi_scientific_name, sort = TRUE)
This example summarizes the source organisms represented in the polymer-entity results. A table like this is often the first step in identifying redundancy, sampling bias, or opportunities to compare orthologous structures across species.
r3d <- asNamespace("r3dmol")
saved_structure <- get_pdb_file(
pdb_id = selected_entry,
filetype = "pdb",
save = TRUE,
path = tempdir(),
verbosity = FALSE
)
r3d$r3dmol() %>%
r3d$m_add_model(data = saved_structure$path, format = "pdb") %>%
r3d$m_set_style(style = r3d$m_style_cartoon(color = "spectrum")) %>%
r3d$m_zoom_to()
This optional chunk demonstrates how rPDBapi fits into
broader R visualization workflows. The package itself focuses on data
access and parsing, while a tool such as r3dmol can be used
to render the retrieved structure in 3D. That separation of
responsibilities is useful because it keeps data access, analysis, and
visualization composable.
The previous sections used text-based and attribute-based search.
rPDBapi also supports sequence, motif, structure, and
chemical searches. These are especially important in structural
bioinformatics, where the biological question is often not “Which entry
mentions a keyword?” but “Which structures resemble this sequence,
motif, fold, or ligand chemistry?”
kinase_motif_sequence <- "VAIKTLKPGTMSPEAFLQEAQVMKKLRHEKLVQLYAVV"
sequence_operator <- SequenceOperator(
sequence = kinase_motif_sequence,
sequence_type = "PROTEIN",
evalue_cutoff = 10,
identity_cutoff = 0.7
)
sequence_operator
#> $evalue_cutoff
#> [1] 10
#>
#> $identity_cutoff
#> [1] 0.7
#>
#> $target
#> [1] "pdb_protein_sequence"
#>
#> $value
#> [1] "VAIKTLKPGTMSPEAFLQEAQVMKKLRHEKLVQLYAVV"
#>
#> attr(,"class")
#> [1] "SequenceOperator" "list"
autoresolve_sequence_type("ATGCGTACGTAGC")
#> [1] "DNA"
autoresolve_sequence_type("AUGCGUACGUAGC")
#> [1] "RNA"
SequenceOperator() constructs a sequence-similarity
search request. This is useful when you have a sequence of interest,
such as a kinase domain fragment, and want to find structurally
characterized homologs in the PDB. The companion function
autoresolve_sequence_type() infers whether a sequence is
DNA, RNA, or protein, which is helpful for interactive workflows and
programmatic pipelines that ingest sequences from external sources.
sequence_hits <- perform_search(
search_operator = sequence_operator,
return_type = "POLYMER_ENTITY",
request_options = RequestOptions(result_start_index = 0, num_results = 5),
verbosity = FALSE
)
sequence_hits
This query asks the RCSB search service for polymer entities similar to the input sequence. The polymer-entity return type is appropriate here because sequence similarity is defined at the entity level rather than at the whole entry level.
prosite_like_motif <- SeqMotifOperator(
pattern = "[LIV]-x-K-[GST]",
sequence_type = "PROTEIN",
pattern_type = "REGEX"
)
prosite_like_motif
#> $value
#> [1] "[LIV]-x-K-[GST]"
#>
#> $pattern_type
#> [1] "regex"
#>
#> $target
#> [1] "pdb_protein_sequence"
#>
#> attr(,"class")
#> [1] "SeqMotifOperator" "list"
SeqMotifOperator() targets local sequence patterns
rather than whole-sequence similarity. This is useful for catalytic
signatures, binding motifs, or short conserved regions that occur across
otherwise diverse proteins.
motif_hits <- perform_search(
search_operator = prosite_like_motif,
return_type = "POLYMER_ENTITY",
request_options = RequestOptions(result_start_index = 0, num_results = 5),
verbosity = FALSE
)
motif_hits
Motif search is especially helpful when a biological question depends on a functional local pattern rather than on full-length homology. In kinase-related work, motifs can help you locate catalytic or regulatory sequence signatures across structural entries.
structure_operator <- StructureOperator(
pdb_entry_id = "4HHB",
assembly_id = 1,
search_mode = "RELAXED_SHAPE_MATCH"
)
structure_operator
#> $value
#> $value$entry_id
#> [1] "4HHB"
#>
#> $value$assembly_id
#> [1] "1"
#>
#>
#> $operator
#> [1] "relaxed_shape_match"
#>
#> attr(,"class")
#> [1] "StructureOperator" "list"
infer_search_service(structure_operator)
#> [1] "structure"
StructureOperator() creates a shape-based structure
search using an existing PDB structure as the template. This is useful
when you want to identify structures with related global geometry, even
if sequence identity is limited. infer_search_service()
confirms that this operator is routed to the structure search backend
rather than to text or full-text search.
structure_hits <- perform_search(
search_operator = structure_operator,
return_type = "ASSEMBLY",
request_options = RequestOptions(result_start_index = 0, num_results = 5),
verbosity = FALSE
)
structure_hits
Returning assemblies is often sensible for shape-based search because biological function frequently depends on the quaternary arrangement of chains rather than only on the deposited asymmetric unit.
atp_like_operator <- ChemicalOperator(
descriptor = "O=P(O)(O)OP(=O)(O)OP(=O)(O)O",
matching_criterion = "fingerprint-similarity"
)
atp_like_operator
#> $value
#> [1] "O=P(O)(O)OP(=O)(O)OP(=O)(O)O"
#>
#> $type
#> [1] "descriptor"
#>
#> $descriptor_type
#> [1] "SMILES"
#>
#> $match_type
#> [1] "fingerprint-similarity"
#>
#> attr(,"class")
#> [1] "ChemicalOperator" "list"
infer_search_service(atp_like_operator)
#> [1] "chemical"
ChemicalOperator() supports ligand-centric workflows by
searching the archive with a SMILES or InChI descriptor. This is useful
when you want to find structures containing chemically similar ligands,
cofactors, or fragments.
chemical_hits <- perform_search(
search_operator = atp_like_operator,
return_type = "CHEMICAL_COMPONENT",
request_options = RequestOptions(result_start_index = 0, num_results = 5),
verbosity = FALSE
)
chemical_hits
Here the return type is CHEMICAL_COMPONENT, which maps
to molecular definitions in the RCSB search API. This lets the search
focus on ligand identities rather than on whole macromolecular
entries.
The operator-based search grammar is one of the most important features of the package. The examples below summarize the text-oriented operators that can be combined in grouped queries.
exact_resolution <- ExactMatchOperator(
attribute = "exptl.method",
value = "X-RAY DIFFRACTION"
)
organism_inclusion <- InOperator(
attribute = "rcsb_entity_source_organism.taxonomy_lineage.name",
value = c("Homo sapiens", "Mus musculus")
)
title_words <- ContainsWordsOperator(
attribute = "struct.title",
value = "protein kinase"
)
title_phrase <- ContainsPhraseOperator(
attribute = "struct.title",
value = "protein kinase"
)
resolution_cutoff <- ComparisonOperator(
attribute = "rcsb_entry_info.resolution_combined",
value = 2.0,
comparison_type = "LESS"
)
resolution_window <- RangeOperator(
attribute = "rcsb_entry_info.resolution_combined",
from_value = 1.0,
to_value = 2.5
)
doi_exists <- ExistsOperator("rcsb_primary_citation.pdbx_database_id_doi")
list(
exact_resolution = exact_resolution,
organism_inclusion = organism_inclusion,
title_words = title_words,
title_phrase = title_phrase,
resolution_cutoff = resolution_cutoff,
resolution_window = resolution_window,
doi_exists = doi_exists
)
#> $exact_resolution
#> $attribute
#> [1] "exptl.method"
#>
#> $value
#> [1] "X-RAY DIFFRACTION"
#>
#> $operator
#> [1] "exact_match"
#>
#> attr(,"class")
#> [1] "ExactMatchOperator" "list"
#>
#> $organism_inclusion
#> $attribute
#> [1] "rcsb_entity_source_organism.taxonomy_lineage.name"
#>
#> $operator
#> [1] "in"
#>
#> $value
#> [1] "Homo sapiens" "Mus musculus"
#>
#> attr(,"class")
#> [1] "InOperator" "list"
#>
#> $title_words
#> $attribute
#> [1] "struct.title"
#>
#> $operator
#> [1] "contains_words"
#>
#> $value
#> [1] "protein kinase"
#>
#> attr(,"class")
#> [1] "ContainsWordsOperator" "list"
#>
#> $title_phrase
#> $attribute
#> [1] "struct.title"
#>
#> $operator
#> [1] "contains_phrase"
#>
#> $value
#> [1] "protein kinase"
#>
#> attr(,"class")
#> [1] "ContainsPhraseOperator" "list"
#>
#> $resolution_cutoff
#> $operator
#> [1] "less"
#>
#> $attribute
#> [1] "rcsb_entry_info.resolution_combined"
#>
#> $value
#> [1] 2
#>
#> attr(,"class")
#> [1] "ComparisonOperator" "list"
#>
#> $resolution_window
#> $operator
#> [1] "range"
#>
#> $attribute
#> [1] "rcsb_entry_info.resolution_combined"
#>
#> $negation
#> [1] FALSE
#>
#> $value
#> $value$from
#> [1] 1
#>
#> $value$to
#> [1] 2.5
#>
#>
#> attr(,"class")
#> [1] "RangeOperator" "list"
#>
#> $doi_exists
#> $attribute
#> [1] "rcsb_primary_citation.pdbx_database_id_doi"
#>
#> $operator
#> [1] "exists"
#>
#> attr(,"class")
#> [1] "ExistsOperator" "list"
These operator constructors do not perform a search by themselves. Instead, they make query intent explicit and composable. That matters when analyses need to be read, reviewed, and reproduced months later.
operator_node <- QueryNode(title_words)
composite_query <- QueryGroup(
queries = list(title_words, resolution_window, doi_exists),
logical_operator = "AND"
)
scored_example <- ScoredResult(entity_id = "4HHB", score = 0.98)
operator_node
#> $type
#> [1] "terminal"
#>
#> $service
#> [1] "text"
#>
#> $parameters
#> $attribute
#> [1] "struct.title"
#>
#> $operator
#> [1] "contains_words"
#>
#> $value
#> [1] "protein kinase"
#>
#> attr(,"class")
#> [1] "ContainsWordsOperator" "list"
composite_query
#> $type
#> [1] "group"
#>
#> $logical_operator
#> [1] "and"
#>
#> $nodes
#> $nodes[[1]]
#> $nodes[[1]]$type
#> [1] "terminal"
#>
#> $nodes[[1]]$service
#> [1] "text"
#>
#> $nodes[[1]]$parameters
#> $attribute
#> [1] "struct.title"
#>
#> $operator
#> [1] "contains_words"
#>
#> $value
#> [1] "protein kinase"
#>
#> attr(,"class")
#> [1] "ContainsWordsOperator" "list"
#>
#>
#> $nodes[[2]]
#> $nodes[[2]]$type
#> [1] "terminal"
#>
#> $nodes[[2]]$service
#> [1] "text"
#>
#> $nodes[[2]]$parameters
#> $operator
#> [1] "range"
#>
#> $attribute
#> [1] "rcsb_entry_info.resolution_combined"
#>
#> $negation
#> [1] FALSE
#>
#> $value
#> $value$from
#> [1] 1
#>
#> $value$to
#> [1] 2.5
#>
#>
#> attr(,"class")
#> [1] "RangeOperator" "list"
#>
#>
#> $nodes[[3]]
#> $nodes[[3]]$type
#> [1] "terminal"
#>
#> $nodes[[3]]$service
#> [1] "text"
#>
#> $nodes[[3]]$parameters
#> $attribute
#> [1] "rcsb_primary_citation.pdbx_database_id_doi"
#>
#> $operator
#> [1] "exists"
#>
#> attr(,"class")
#> [1] "ExistsOperator" "list"
scored_example
#> $entity_id
#> [1] "4HHB"
#>
#> $score
#> [1] 0.98
QueryNode() and QueryGroup() are the glue
that turn independent operator objects into a full search graph.
ScoredResult() is a small utility that represents the shape
of a scored hit and is useful for result handling or for teaching the
output model used by structure-search APIs.
scored_structure_hits <- perform_search(
search_operator = structure_operator,
return_type = "ASSEMBLY",
request_options = RequestOptions(result_start_index = 0, num_results = 3),
return_with_scores = TRUE,
verbosity = FALSE
)
scored_structure_hits
class(scored_structure_hits)
This example shows the difference between returning plain identifiers and returning scored hits. In similarity-oriented workflows, the score itself can be analytically useful because it helps rank follow-up candidates before any metadata are fetched. A common pattern is to inspect scored results first, choose a cutoff, and only then retrieve metadata for the retained identifiers.
# Pattern: build small reusable operators first
title_filter <- ContainsPhraseOperator("struct.title", "protein kinase")
resolution_filter <- ComparisonOperator(
"rcsb_entry_info.resolution_combined",
2.5,
"LESS_OR_EQUAL"
)
# Wrap or combine them only when the biological question is clear
query_graph <- QueryGroup(
queries = list(
QueryNode(title_filter),
QueryNode(resolution_filter)
),
logical_operator = "AND"
)
This is the main reason to treat operator construction as a separate step. A query graph can be assembled gradually, reviewed independently of the network request, and reused across multiple searches or result types.
query_search() is intentionally simpler than
perform_search(), but it still supports several specialized
query types as well as low-level request overrides.
# PubMed-linked structures
query_search(search_term = 27499440, query_type = "PubmedIdQuery")
# Organism/taxonomy search
query_search(search_term = "9606", query_type = "TreeEntityQuery")
# Experimental method search
query_search(search_term = "X-RAY DIFFRACTION", query_type = "ExpTypeQuery")
# Author search
query_search(search_term = "Kuriyan, J.", query_type = "AdvancedAuthorQuery")
# UniProt-linked entries
query_search(search_term = "P31749", query_type = "uniprot")
# PFAM-linked entries
query_search(search_term = "PF00069", query_type = "pfam")
These convenience modes are useful when the search criterion maps directly to a common biological identifier or curation field. They are less flexible than a fully operator-based query, but faster to write for routine tasks.
custom_scan_params <- list(
request_options = list(
paginate = list(start = 0, rows = 5),
return_all_hits = FALSE
)
)
custom_scan_params
#> $request_options
#> $request_options$paginate
#> $request_options$paginate$start
#> [1] 0
#>
#> $request_options$paginate$rows
#> [1] 5
#>
#>
#> $request_options$return_all_hits
#> [1] FALSE
scan_params lets you override the request body that
query_search() sends to the search API. This is useful when
you want lightweight access to custom pagination or request options
without switching fully to perform_search().
limited_kinase_hits <- query_search(
search_term = "protein kinase",
scan_params = custom_scan_params
)
limited_kinase_hits
This example illustrates the practical use of
scan_params: constrain the result set while still using the
simpler query helper.
The data_fetcher() interface supports more than entry
and polymer-entity data. The main supported data types are:
ENTRYASSEMBLYPOLYMER_ENTITYBRANCHED_ENTITYNONPOLYMER_ENTITYPOLYMER_ENTITY_INSTANCEBRANCHED_ENTITY_INSTANCENONPOLYMER_ENTITY_INSTANCECHEMICAL_COMPONENTThis breadth matters because structural records are hierarchical. Different questions belong to different levels: entry-level methods, assembly-level symmetry, entity-level taxonomy, instance-level chain annotations, and component-level ligand chemistry.
base_properties <- list(
rcsb_entry_info = c("resolution_combined"),
exptl = c("method")
)
extended_properties <- add_property(list(
rcsb_entry_info = c("molecular_weight", "resolution_combined"),
struct = c("title")
))
base_properties
#> $rcsb_entry_info
#> [1] "resolution_combined"
#>
#> $exptl
#> [1] "method"
extended_properties
#> $rcsb_entry_info
#> [1] "molecular_weight" "resolution_combined"
#>
#> $struct
#> [1] "title"
add_property() helps construct or merge property
specifications without duplicating subfields. This is especially useful
in interactive analyses, where you may start with a minimal query and
then progressively request additional annotations.
property_workflow <- add_property(list(
rcsb_id = list(),
struct = c("title"),
rcsb_entry_info = c("resolution_combined")
))
property_workflow <- add_property(list(
rcsb_entry_info = c("molecular_weight", "resolution_combined"),
exptl = c("method")
))
property_workflow
#> $rcsb_entry_info
#> [1] "molecular_weight" "resolution_combined"
#>
#> $exptl
#> [1] "method"
validate_properties(property_workflow, data_type = "ENTRY", strict = FALSE)
#> $unknown_fields
#> character(0)
#>
#> $unknown_subfields
#> list()
This pattern is useful because GraphQL property lists tend to grow as an analysis becomes more specific. Building them incrementally makes it easier to keep a compact initial query, add only the fields that become necessary, and check that the evolving specification still matches the expected schema.
ligand_properties <- list(
rcsb_id = list(),
chem_comp = c("id", "name", "formula", "formula_weight", "type"),
rcsb_chem_comp_info = c("initial_release_date")
)
ligand_properties
#> $rcsb_id
#> list()
#>
#> $chem_comp
#> [1] "id" "name" "formula" "formula_weight"
#> [5] "type"
#>
#> $rcsb_chem_comp_info
#> [1] "initial_release_date"
This property specification targets chemical components rather than whole structures. That distinction is important when the biological focus is on ligands, cofactors, inhibitors, or bound metabolites.
chemical_component_df <- data_fetcher(
id = head(chemical_hits, 3),
data_type = "CHEMICAL_COMPONENT",
properties = ligand_properties,
return_as_dataframe = TRUE
)
chemical_component_df
The resulting table can be used to compare ligand formulas, molecular weights, and release histories. This is often useful in medicinal chemistry and structure-based design workflows.
ligand_object <- as_rpdb_chemical_component(
chemical_component_df,
metadata = list(query = "ATP-like chemical components")
)
extract_ligand_table(ligand_object)
extract_ligand_table() keeps the most analysis-relevant
chemical-component columns in a compact form. That is useful when ligand
retrieval is only one part of a broader workflow and you want a small,
stable table for downstream joins or ranking.
atp_description <- describe_chemical("ATP")
dplyr::tibble(
chem_id = "ATP",
name = purrr::pluck(atp_description, "chem_comp", "name", .default = NA_character_),
formula = purrr::pluck(atp_description, "chem_comp", "formula", .default = NA_character_),
formula_weight = purrr::pluck(atp_description, "chem_comp", "formula_weight", .default = NA),
smiles = purrr::pluck(atp_description, "rcsb_chem_comp_descriptor", "smiles", .default = NA_character_)
)
describe_chemical() provides a direct route to detailed
ligand information for a single chemical component. It complements
data_fetcher() by supporting a focused, ligand-centric
lookup.
# Polymer chain instances
data_fetcher(
id = c("4HHB.A"),
data_type = "POLYMER_ENTITY_INSTANCE",
properties = list(rcsb_id = list())
)
# Branched entity instances
data_fetcher(
id = c("1ABC_2"),
data_type = "BRANCHED_ENTITY_INSTANCE",
properties = list(rcsb_id = list())
)
# Non-polymer entity instances
data_fetcher(
id = c("3PQR_5"),
data_type = "NONPOLYMER_ENTITY_INSTANCE",
properties = list(rcsb_id = list())
)
Instance-level retrieval is relevant when chain-level or site-specific annotations matter. The exact identifier format depends on the corresponding RCSB data type and record level.
The package exposes lower-level functions for users who need full control over HTTP requests, URLs, or response parsing.
entry_url <- get_pdb_api_url("core/entry/", "4HHB")
chem_url <- get_pdb_api_url("core/chemcomp/", "ATP")
entry_url
#> [1] "https://data.rcsb.org/rest/v1/core/entry/4HHB"
chem_url
#> [1] "https://data.rcsb.org/rest/v1/core/chemcomp/ATP"
get_pdb_api_url() constructs endpoint-specific URLs.
This is a small utility, but it makes low-level workflows clearer and
reduces hard-coded URL strings in custom scripts.
# Manual request lifecycle
url <- get_pdb_api_url("core/entry/", "4HHB")
response <- send_api_request(url, verbosity = FALSE)
handle_api_errors(response, url)
payload <- parse_response(response, format = "json")
This low-level lifecycle is useful when you are developing a new helper, debugging an endpoint transition, or comparing package behavior with the raw RCSB API. It also makes clear where URL construction, HTTP transport, status checking, and parsing are separated inside the package.
entry_response <- send_api_request(entry_url, verbosity = FALSE)
handle_api_errors(entry_response, entry_url)
entry_payload <- parse_response(entry_response, format = "json")
names(entry_payload)[1:5]
These functions expose the package’s low-level request stack.
send_api_request() handles the HTTP request,
handle_api_errors() checks the returned status, and
parse_response() converts the body into an R object. This
layer is useful when you need to debug endpoint behavior or prototype a
new helper around the RCSB REST API.
mini_graphql <- generate_json_query(
ids = kinase_entry_ids[1:2],
data_type = "ENTRY",
properties = list(rcsb_id = list(), struct = c("title"))
)
mini_graphql_response <- search_graphql(mini_graphql)
str(mini_graphql_response, max.level = 2)
search_graphql() is the low-level GraphQL entry point.
It is useful when you want to inspect the raw content returned by the
RCSB GraphQL service before it is normalized by
fetch_data() or flattened by
return_data_as_dataframe().
One of the notable features of the current package version is that core functions return typed objects. This improves programmatic safety because code can distinguish search identifiers, scored results, raw responses, and flattened data frames.
list(
query_search_class = class(query_search("kinase")),
perform_search_class = class(
perform_search(DefaultOperator("kinase"), verbosity = FALSE)
),
perform_search_scores_class = class(
perform_search(
DefaultOperator("kinase"),
return_with_scores = TRUE,
verbosity = FALSE
)
)
)
The classes shown above make return semantics explicit. In a larger analysis pipeline, this reduces ambiguity and makes it easier to validate assumptions at each stage of data retrieval.
raw_entry_response <- data_fetcher(
id = kinase_entry_ids[1:2],
data_type = "ENTRY",
properties = list(rcsb_id = list()),
return_as_dataframe = FALSE
)
tidy_entry_response <- data_fetcher(
id = kinase_entry_ids[1:2],
data_type = "ENTRY",
properties = list(rcsb_id = list()),
return_as_dataframe = TRUE
)
class(raw_entry_response)
class(tidy_entry_response)
This example emphasizes the dual output model of
data_fetcher(): retain the nested payload when structure
matters, or request a data frame when analysis and integration matter
more.
list(
entry_object_class = class(as_rpdb_entry(kinase_metadata)),
assembly_object_class = class(as_rpdb_assembly(kinase_assemblies)),
polymer_object_class = class(as_rpdb_polymer_entity(kinase_polymer_metadata)),
structure_object_class = class(as_rpdb_structure(kinase_structure)),
batch_provenance_names = names(attr(kinase_batch, "provenance"))
)
These richer object wrappers are deliberately lightweight. They preserve the underlying data while attaching a semantically meaningful class, which makes it easier to define helper methods and to branch on object type in larger structural analysis pipelines.
local_entry_object <- as_rpdb_entry(
data.frame(
rcsb_id = "4HHB",
method = "X-RAY DIFFRACTION",
resolution_combined = "1.74",
stringsAsFactors = FALSE
),
metadata = list(source = "local method demo")
)
print(local_entry_object)
#> <rPDBapi_entry> with data class: data.frame
dplyr::as_tibble(local_entry_object)
#> # A tibble: 1 × 3
#> rcsb_id method resolution_combined
#> <chr> <chr> <chr>
#> 1 4HHB X-RAY DIFFRACTION 1.74
This example makes the object behavior explicit. The custom print
method gives a concise summary of the wrapped object, while the
as_tibble() method provides an immediate path back to a
standard tabular workflow. That combination is the main point of the
object model: preserve semantic type information without making
downstream manipulation cumbersome.
invalid_property_result <- tryCatch(
validate_properties(
properties = list(unknown_field = c("x")),
data_type = "ENTRY",
strict = TRUE
),
rPDBapi_error_invalid_input = function(e) e
)
invalid_fetch_result <- tryCatch(
data_fetcher(
id = character(0),
data_type = "ENTRY",
properties = list(rcsb_id = list())
),
rPDBapi_error_invalid_input = function(e) e
)
list(
invalid_property_class = class(invalid_property_result),
invalid_property_message = conditionMessage(invalid_property_result),
invalid_fetch_class = class(invalid_fetch_result),
invalid_fetch_message = conditionMessage(invalid_fetch_result)
)
#> $invalid_property_class
#> [1] "rPDBapi_error_invalid_input" "rPDBapi_error"
#> [3] "error" "condition"
#>
#> $invalid_property_message
#> [1] "Unknown properties for data_type 'ENTRY': unknown_field"
#>
#> $invalid_fetch_class
#> [1] "rPDBapi_error_invalid_input" "rPDBapi_error"
#> [3] "error" "condition"
#>
#> $invalid_fetch_message
#> [1] "Invalid input: 'id' must not be NULL or empty."
These examples show how typed errors support defensive programming.
Instead of matching raw error text, a calling workflow can branch on the
condition class and decide whether to stop, retry, skip a record, or log
the problem for later review. That is particularly valuable when
rPDBapi is used inside larger automated pipelines.
The table below maps every exported function to its primary role in the package. This is not a replacement for the individual help pages, but it does make the full surface area of the package explicit inside a single tutorial document.
export_reference <- data.frame(
Function = c(
"query_search", "perform_search", "DefaultOperator", "ExactMatchOperator",
"InOperator", "ContainsWordsOperator", "ContainsPhraseOperator",
"ComparisonOperator", "RangeOperator", "ExistsOperator",
"SequenceOperator", "autoresolve_sequence_type", "SeqMotifOperator",
"StructureOperator", "ChemicalOperator", "QueryNode", "QueryGroup",
"RequestOptions", "ScoredResult", "infer_search_service",
"infer_id_type", "parse_rcsb_id", "build_entry_id", "build_assembly_id",
"build_entity_id", "build_instance_id", "add_property",
"list_rcsb_fields", "search_rcsb_fields", "validate_properties",
"generate_json_query", "search_graphql", "fetch_data",
"return_data_as_dataframe", "data_fetcher", "data_fetcher_batch",
"cache_info", "clear_rpdbapi_cache", "get_info", "find_results",
"find_papers", "describe_chemical", "get_fasta_from_rcsb_entry",
"get_pdb_file", "get_pdb_api_url", "send_api_request",
"handle_api_errors", "parse_response", "as_rpdb_entry",
"as_rpdb_assembly", "as_rpdb_polymer_entity",
"as_rpdb_chemical_component", "as_rpdb_structure",
"summarize_entries", "summarize_assemblies",
"extract_taxonomy_table", "extract_ligand_table",
"extract_calpha_coordinates", "join_structure_sequence"
),
Role = c(
"High-level convenience search helper",
"Operator-based search engine",
"Full-text search operator",
"Exact attribute match operator",
"Set-membership operator",
"Word containment operator",
"Phrase containment operator",
"Numeric/date comparison operator",
"Range filter operator",
"Attribute existence operator",
"Sequence similarity search operator",
"Automatic DNA/RNA/protein detection",
"Sequence motif search operator",
"Structure similarity search operator",
"Chemical descriptor search operator",
"Wrap one operator as a query node",
"Combine nodes with AND/OR logic",
"Pagination and sorting controls",
"Represent a scored hit",
"Infer backend service from operator",
"Infer identifier level from an ID string",
"Parse an identifier into structured components",
"Normalize or build entry identifiers",
"Build assembly identifiers",
"Build entity identifiers",
"Build instance or chain identifiers",
"Merge/extend GraphQL property lists",
"List known retrievable fields by data type",
"Search the built-in field registry",
"Validate a property list against the field registry",
"Build a GraphQL query string",
"Low-level GraphQL request helper",
"Normalize validated GraphQL payloads",
"Flatten nested payloads into data frames",
"High-level metadata fetcher",
"Batch metadata fetcher with retry and provenance",
"Inspect batch-cache contents",
"Clear on-disk cache entries",
"Retrieve full entry metadata",
"Extract one field across search hits",
"Extract primary citation titles",
"Retrieve ligand/chemical-component details",
"Retrieve FASTA sequences",
"Download and parse structure files",
"Build REST endpoint URLs",
"Send low-level GET/POST requests",
"Check HTTP status and stop on error",
"Parse JSON or text responses",
"Wrap entry data in a typed object",
"Wrap assembly data in a typed object",
"Wrap polymer-entity data in a typed object",
"Wrap chemical-component data in a typed object",
"Wrap structure data in a typed object",
"Summarize entry-level metadata",
"Summarize assembly-level metadata",
"Extract taxonomy-focused columns",
"Extract ligand-focused columns",
"Extract C-alpha coordinates",
"Join sequence summaries to chain coordinates"
),
stringsAsFactors = FALSE
)
knitr::kable(export_reference, align = c("l", "l"))
| Function | Role |
|---|---|
| query_search | High-level convenience search helper |
| perform_search | Operator-based search engine |
| DefaultOperator | Full-text search operator |
| ExactMatchOperator | Exact attribute match operator |
| InOperator | Set-membership operator |
| ContainsWordsOperator | Word containment operator |
| ContainsPhraseOperator | Phrase containment operator |
| ComparisonOperator | Numeric/date comparison operator |
| RangeOperator | Range filter operator |
| ExistsOperator | Attribute existence operator |
| SequenceOperator | Sequence similarity search operator |
| autoresolve_sequence_type | Automatic DNA/RNA/protein detection |
| SeqMotifOperator | Sequence motif search operator |
| StructureOperator | Structure similarity search operator |
| ChemicalOperator | Chemical descriptor search operator |
| QueryNode | Wrap one operator as a query node |
| QueryGroup | Combine nodes with AND/OR logic |
| RequestOptions | Pagination and sorting controls |
| ScoredResult | Represent a scored hit |
| infer_search_service | Infer backend service from operator |
| infer_id_type | Infer identifier level from an ID string |
| parse_rcsb_id | Parse an identifier into structured components |
| build_entry_id | Normalize or build entry identifiers |
| build_assembly_id | Build assembly identifiers |
| build_entity_id | Build entity identifiers |
| build_instance_id | Build instance or chain identifiers |
| add_property | Merge/extend GraphQL property lists |
| list_rcsb_fields | List known retrievable fields by data type |
| search_rcsb_fields | Search the built-in field registry |
| validate_properties | Validate a property list against the field registry |
| generate_json_query | Build a GraphQL query string |
| search_graphql | Low-level GraphQL request helper |
| fetch_data | Normalize validated GraphQL payloads |
| return_data_as_dataframe | Flatten nested payloads into data frames |
| data_fetcher | High-level metadata fetcher |
| data_fetcher_batch | Batch metadata fetcher with retry and provenance |
| cache_info | Inspect batch-cache contents |
| clear_rpdbapi_cache | Clear on-disk cache entries |
| get_info | Retrieve full entry metadata |
| find_results | Extract one field across search hits |
| find_papers | Extract primary citation titles |
| describe_chemical | Retrieve ligand/chemical-component details |
| get_fasta_from_rcsb_entry | Retrieve FASTA sequences |
| get_pdb_file | Download and parse structure files |
| get_pdb_api_url | Build REST endpoint URLs |
| send_api_request | Send low-level GET/POST requests |
| handle_api_errors | Check HTTP status and stop on error |
| parse_response | Parse JSON or text responses |
| as_rpdb_entry | Wrap entry data in a typed object |
| as_rpdb_assembly | Wrap assembly data in a typed object |
| as_rpdb_polymer_entity | Wrap polymer-entity data in a typed object |
| as_rpdb_chemical_component | Wrap chemical-component data in a typed object |
| as_rpdb_structure | Wrap structure data in a typed object |
| summarize_entries | Summarize entry-level metadata |
| summarize_assemblies | Summarize assembly-level metadata |
| extract_taxonomy_table | Extract taxonomy-focused columns |
| extract_ligand_table | Extract ligand-focused columns |
| extract_calpha_coordinates | Extract C-alpha coordinates |
| join_structure_sequence | Join sequence summaries to chain coordinates |
This table is intended as a package navigation aid. It makes it easier to identify whether a task belongs to searching, retrieval, parsing, or lower-level API control before you start writing a workflow.
The next block gives a compact usage sketch for every exported function. These examples are deliberately short and grouped by role so that users can quickly find a starting pattern.
# Search helpers
query_search("kinase")
perform_search(DefaultOperator("kinase"), verbosity = FALSE)
# Text and attribute operators
DefaultOperator("kinase")
ExactMatchOperator("exptl.method", "X-RAY DIFFRACTION")
InOperator("rcsb_entity_source_organism.taxonomy_lineage.name", c("Homo sapiens", "Mus musculus"))
ContainsWordsOperator("struct.title", "protein kinase")
ContainsPhraseOperator("struct.title", "protein kinase")
ComparisonOperator("rcsb_entry_info.resolution_combined", 2.0, "LESS")
RangeOperator("rcsb_entry_info.resolution_combined", 1.0, 2.5)
ExistsOperator("rcsb_primary_citation.pdbx_database_id_doi")
# Specialized operators
SequenceOperator("MVLSPADKTNVKAAW", sequence_type = "PROTEIN")
autoresolve_sequence_type("ATGCGTACGTAGC")
SeqMotifOperator("[LIV]-x-K-[GST]", "PROTEIN", "REGEX")
StructureOperator("4HHB", assembly_id = 1, search_mode = "RELAXED_SHAPE_MATCH")
ChemicalOperator("C1=CC=CC=C1", matching_criterion = "graph-strict")
# Query composition
QueryNode(DefaultOperator("kinase"))
QueryGroup(list(DefaultOperator("kinase"), ExistsOperator("rcsb_primary_citation.title")), "AND")
RequestOptions(result_start_index = 0, num_results = 10)
ScoredResult("4HHB", 0.98)
infer_search_service(StructureOperator("4HHB"))
infer_id_type(c("4HHB", "4HHB-1", "4HHB_1", "4HHB.A", "ATP"))
parse_rcsb_id("4HHB-1")
build_entry_id("4HHB")
build_assembly_id("4HHB", 1)
build_entity_id("4HHB", 1)
build_instance_id("4HHB", "A")
# Metadata helpers
add_property(list(rcsb_entry_info = c("resolution_combined")))
list_rcsb_fields("ENTRY")
search_rcsb_fields("resolution", data_type = "ENTRY")
validate_properties(
list(rcsb_id = list(), rcsb_entry_info = c("resolution_combined")),
data_type = "ENTRY",
strict = TRUE
)
generate_json_query(c("4HHB"), "ENTRY", list(rcsb_id = list(), struct = c("title")))
search_graphql(generate_json_query(c("4HHB"), "ENTRY", list(rcsb_id = list())))
fetch_data(generate_json_query(c("4HHB"), "ENTRY", list(rcsb_id = list())), "ENTRY", "4HHB")
return_data_as_dataframe(
fetch_data(generate_json_query(c("4HHB"), "ENTRY", list(rcsb_id = list())), "ENTRY", "4HHB"),
"ENTRY",
"4HHB"
)
data_fetcher("4HHB", "ENTRY", list(rcsb_id = list(), struct = c("title")))
data_fetcher_batch(
c("4HHB", "1CRN"),
"ENTRY",
list(rcsb_id = list(), struct = c("title")),
batch_size = 1,
cache = FALSE
)
cache_info()
clear_rpdbapi_cache()
get_info("4HHB")
find_results("kinase", field = "struct_keywords")
find_papers("kinase", max_results = 3)
describe_chemical("ATP")
get_fasta_from_rcsb_entry("4HHB")
# Files and low-level HTTP
get_pdb_file("4HHB", filetype = "cif", verbosity = FALSE)
get_pdb_api_url("core/entry/", "4HHB")
resp <- send_api_request(get_pdb_api_url("core/entry/", "4HHB"), verbosity = FALSE)
handle_api_errors(resp, get_pdb_api_url("core/entry/", "4HHB"))
parse_response(resp, format = "json")
# Object wrappers and analysis helpers
as_rpdb_entry(data.frame(rcsb_id = "4HHB"))
as_rpdb_assembly(data.frame(rcsb_id = "4HHB-1"))
as_rpdb_polymer_entity(data.frame(rcsb_id = "4HHB_1"))
as_rpdb_chemical_component(data.frame(rcsb_id = "ATP"))
as_rpdb_structure(get_pdb_file("4HHB", filetype = "cif", verbosity = FALSE))
summarize_entries(data.frame(method = "X-RAY DIFFRACTION", resolution_combined = "1.8"))
summarize_assemblies(data.frame(oligomeric_count = "2", symbol = "C2"))
extract_taxonomy_table(data.frame(rcsb_id = "4HHB_1", ncbi_taxonomy_id = "9606"))
extract_ligand_table(data.frame(rcsb_id = "ATP", formula_weight = "507.18"))
extract_calpha_coordinates(get_pdb_file("4HHB", filetype = "cif", verbosity = FALSE))
join_structure_sequence(
get_pdb_file("4HHB", filetype = "cif", verbosity = FALSE),
get_fasta_from_rcsb_entry("4HHB")
)
This appendix is intentionally compact. Its purpose is not to replace the narrative examples above, but to ensure that every exported function has an immediately visible calling pattern in the vignette.
One source of confusion in the RCSB ecosystem is that different
endpoints expect different identifier types. The table below summarizes
the levels supported by data_fetcher() and the search
return types used in perform_search().
id_reference <- data.frame(
Data_or_Return_Type = c(
"ENTRY", "ASSEMBLY", "POLYMER_ENTITY", "BRANCHED_ENTITY",
"NONPOLYMER_ENTITY", "POLYMER_ENTITY_INSTANCE",
"BRANCHED_ENTITY_INSTANCE", "NONPOLYMER_ENTITY_INSTANCE",
"CHEMICAL_COMPONENT"
),
Typical_ID_Format = c(
"4-character PDB ID, e.g. 4HHB",
"Entry plus assembly ID, e.g. 4HHB-1",
"Entry plus entity ID, e.g. 4HHB_1",
"Entry plus branched entity ID",
"Entry plus nonpolymer entity ID, e.g. 3PQR_5",
"Instance or chain-level identifier, endpoint-specific",
"Instance-level identifier, endpoint-specific",
"Instance-level identifier, endpoint-specific",
"Chemical component ID, e.g. ATP"
),
Typical_Use = c(
"Whole-structure metadata",
"Biological assembly and symmetry",
"Entity-level taxonomy or sequence annotations",
"Glycan/branched entity records",
"Ligand records within structures",
"Chain-specific annotations",
"Branched entity instance records",
"Ligand instance records",
"Ligand chemistry and descriptors"
),
stringsAsFactors = FALSE
)
knitr::kable(id_reference, align = c("l", "l", "l"))
| Data_or_Return_Type | Typical_ID_Format | Typical_Use |
|---|---|---|
| ENTRY | 4-character PDB ID, e.g. 4HHB | Whole-structure metadata |
| ASSEMBLY | Entry plus assembly ID, e.g. 4HHB-1 | Biological assembly and symmetry |
| POLYMER_ENTITY | Entry plus entity ID, e.g. 4HHB_1 | Entity-level taxonomy or sequence annotations |
| BRANCHED_ENTITY | Entry plus branched entity ID | Glycan/branched entity records |
| NONPOLYMER_ENTITY | Entry plus nonpolymer entity ID, e.g. 3PQR_5 | Ligand records within structures |
| POLYMER_ENTITY_INSTANCE | Instance or chain-level identifier, endpoint-specific | Chain-specific annotations |
| BRANCHED_ENTITY_INSTANCE | Instance-level identifier, endpoint-specific | Branched entity instance records |
| NONPOLYMER_ENTITY_INSTANCE | Instance-level identifier, endpoint-specific | Ligand instance records |
| CHEMICAL_COMPONENT | Chemical component ID, e.g. ATP | Ligand chemistry and descriptors |
The precise identifier syntax for instance-level records depends on
the RCSB schema and endpoint, but the key conceptual point is that the
package supports multiple biological levels and expects identifiers that
match those levels. The identifier helpers introduced in the current
version make this easier to manage explicitly:
infer_id_type() classifies common patterns,
parse_rcsb_id() decomposes them, and the
build_*_id() functions generate normalized identifiers
programmatically.
The package uses return classes as lightweight contracts. The most important ones are summarized here.
contract_reference <- data.frame(
Function = c(
"query_search(return_type = 'entry')",
"query_search(other return_type)",
"perform_search()",
"perform_search(return_with_scores = TRUE)",
"perform_search(return_raw_json_dict = TRUE)",
"fetch_data()",
"data_fetcher_batch(return_as_dataframe = TRUE)",
"data_fetcher(return_as_dataframe = TRUE)",
"data_fetcher(return_as_dataframe = FALSE)",
"as_rpdb_entry()",
"as_rpdb_assembly()",
"as_rpdb_polymer_entity()",
"as_rpdb_chemical_component()",
"as_rpdb_structure()"
),
Return_Class = c(
"rPDBapi_query_ids",
"rPDBapi_query_response",
"rPDBapi_search_ids",
"rPDBapi_search_scores",
"rPDBapi_search_raw_response",
"rPDBapi_fetch_response",
"rPDBapi_dataframe",
"rPDBapi_dataframe",
"rPDBapi_fetch_response",
"rPDBapi_entry",
"rPDBapi_assembly",
"rPDBapi_polymer_entity",
"rPDBapi_chemical_component",
"rPDBapi_structure"
),
Meaning = c(
"Identifier vector from query_search()",
"Parsed query_search payload",
"Identifier vector from perform_search()",
"Scored search results",
"Raw JSON-like search payload",
"Validated GraphQL fetch payload",
"Flattened batch result with provenance metadata",
"Flattened analysis-ready table",
"Nested validated fetch payload",
"Typed entry wrapper around retrieved data",
"Typed assembly wrapper around retrieved data",
"Typed polymer-entity wrapper around retrieved data",
"Typed chemical-component wrapper around retrieved data",
"Typed structure wrapper around retrieved data"
),
stringsAsFactors = FALSE
)
knitr::kable(contract_reference, align = c("l", "l", "l"))
| Function | Return_Class | Meaning |
|---|---|---|
| query_search(return_type = ‘entry’) | rPDBapi_query_ids | Identifier vector from query_search() |
| query_search(other return_type) | rPDBapi_query_response | Parsed query_search payload |
| perform_search() | rPDBapi_search_ids | Identifier vector from perform_search() |
| perform_search(return_with_scores = TRUE) | rPDBapi_search_scores | Scored search results |
| perform_search(return_raw_json_dict = TRUE) | rPDBapi_search_raw_response | Raw JSON-like search payload |
| fetch_data() | rPDBapi_fetch_response | Validated GraphQL fetch payload |
| data_fetcher_batch(return_as_dataframe = TRUE) | rPDBapi_dataframe | Flattened batch result with provenance metadata |
| data_fetcher(return_as_dataframe = TRUE) | rPDBapi_dataframe | Flattened analysis-ready table |
| data_fetcher(return_as_dataframe = FALSE) | rPDBapi_fetch_response | Nested validated fetch payload |
| as_rpdb_entry() | rPDBapi_entry | Typed entry wrapper around retrieved data |
| as_rpdb_assembly() | rPDBapi_assembly | Typed assembly wrapper around retrieved data |
| as_rpdb_polymer_entity() | rPDBapi_polymer_entity | Typed polymer-entity wrapper around retrieved data |
| as_rpdb_chemical_component() | rPDBapi_chemical_component | Typed chemical-component wrapper around retrieved data |
| as_rpdb_structure() | rPDBapi_structure | Typed structure wrapper around retrieved data |
These classes are useful when writing wrappers, tests, or pipelines that need to branch on the kind of object returned by the package.
The package also uses typed errors in several important places. Users do not need to memorize these classes for normal interactive work, but they are useful for robust scripting and package development.
error_guidance <- data.frame(
Scenario = c(
"Malformed search response",
"Unsupported return-type mapping",
"Invalid input to search/fetch helper",
"Unknown property or subproperty in strict mode",
"Batch retrieval failure after retries",
"HTTP failure",
"Response parsing failure"
),
Typical_Class_or_Source = c(
"rPDBapi_error_malformed_response",
"rPDBapi_error_unsupported_mapping",
"rPDBapi_error_invalid_input",
"validate_properties() / generate_json_query()",
"data_fetcher_batch()",
"handle_api_errors() / send_api_request()",
"parse_response()"
),
stringsAsFactors = FALSE
)
knitr::kable(error_guidance, align = c("l", "l"))
| Scenario | Typical_Class_or_Source |
|---|---|
| Malformed search response | rPDBapi_error_malformed_response |
| Unsupported return-type mapping | rPDBapi_error_unsupported_mapping |
| Invalid input to search/fetch helper | rPDBapi_error_invalid_input |
| Unknown property or subproperty in strict mode | validate_properties() / generate_json_query() |
| Batch retrieval failure after retries | data_fetcher_batch() |
| HTTP failure | handle_api_errors() / send_api_request() |
| Response parsing failure | parse_response() |
In practice, these classes matter when you want to distinguish network failures from schema mismatches or user-input errors. That distinction is particularly important in automated structural bioinformatics pipelines that may run over many identifiers.
Programmatic structure retrieval is most useful when the search logic and the retrieved identifiers are stored alongside the analysis. A practical workflow is to save:
analysis_manifest <- list(
live_examples = live_examples,
package_version = as.character(utils::packageVersion("rPDBapi")),
query = kinase_query,
requested_entry_fields = entry_properties,
strict_property_validation = getOption("rPDBapi.strict_property_validation", FALSE),
built_ids = list(
entry = build_entry_id("4HHB"),
assembly = build_assembly_id("4HHB", 1),
entity = build_entity_id("4HHB", 1),
instance = build_instance_id("4HHB", "A")
),
batch_provenance_example = if (live_examples) {
attr(kinase_batch, "provenance")
} else {
NULL
}
)
str(analysis_manifest, max.level = 2)
#> List of 7
#> $ live_examples : logi FALSE
#> $ package_version : chr "3.0.0"
#> $ query :List of 3
#> ..$ type : chr "group"
#> ..$ logical_operator: chr "and"
#> ..$ nodes :List of 3
#> $ requested_entry_fields :List of 6
#> ..$ rcsb_id : list()
#> ..$ struct : chr "title"
#> ..$ struct_keywords : chr "pdbx_keywords"
#> ..$ exptl : chr "method"
#> ..$ rcsb_entry_info : chr [1:2] "molecular_weight" "resolution_combined"
#> ..$ rcsb_accession_info: chr "initial_release_date"
#> $ strict_property_validation: logi FALSE
#> $ built_ids :List of 4
#> ..$ entry : chr "4HHB"
#> ..$ assembly: chr "4HHB-1"
#> ..$ entity : chr "4HHB_1"
#> ..$ instance: chr "4HHB.A"
#> $ batch_provenance_example : NULL
This manifest is a simple example of how to preserve the logic of an
analysis. Because the search operators and requested fields are explicit
R objects, they can be saved with saveRDS() and reused
later. That is a better long-term strategy than relying on manual notes
about which website filters were used. When batch retrieval is part of
the workflow, the provenance attribute from
data_fetcher_batch() provides an additional audit trail for
how the data were obtained.
rPDBapi supports an end-to-end workflow for structural
bioinformatics in R: search the archive, refine the result set with
explicit operators, validate properties against known schema fields,
work across identifier levels, retrieve entry-, entity-, and
assembly-level metadata, scale retrieval with batch and cache-aware
helpers, convert nested responses into tidy data frames or typed
objects, download coordinate files, and integrate the results with
analysis and visualization packages. This workflow is useful not only
for exploratory access to the PDB, but also for reproducible, scriptable
analyses that can be revised and rerun as biological questions
evolve.
sessionInfo()
#> R version 4.5.2 (2025-10-31)
#> Platform: aarch64-apple-darwin20
#> Running under: macOS Sequoia 15.7.3
#>
#> Matrix products: default
#> BLAS: /System/Library/Frameworks/Accelerate.framework/Versions/A/Frameworks/vecLib.framework/Versions/A/libBLAS.dylib
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.1
#>
#> locale:
#> [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
#>
#> time zone: Europe/Istanbul
#> tzcode source: internal
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] dplyr_1.2.0 rPDBapi_3.0.0
#>
#> loaded via a namespace (and not attached):
#> [1] vctrs_0.7.1 httr_1.4.8 cli_3.6.5 knitr_1.51
#> [5] rlang_1.1.7 xfun_0.56 otel_0.2.0 purrr_1.2.1
#> [9] bio3d_2.4-5 generics_0.1.4 jsonlite_2.0.0 glue_1.8.0
#> [13] htmltools_0.5.9 sass_0.4.10 rmarkdown_2.30 grid_4.5.2
#> [17] tibble_3.3.1 evaluate_1.0.5 jquerylib_0.1.4 fastmap_1.2.0
#> [21] yaml_2.3.12 lifecycle_1.0.5 compiler_4.5.2 pkgconfig_2.0.3
#> [25] Rcpp_1.1.1 rstudioapi_0.18.0 digest_0.6.39 R6_2.6.1
#> [29] utf8_1.2.6 tidyselect_1.2.1 pillar_1.11.1 parallel_4.5.2
#> [33] magrittr_2.0.4 bslib_0.10.0 tools_4.5.2 xml2_1.5.2
#> [37] cachem_1.1.0
Session information is included so that rendered results can be tied to a specific R and package environment. This is a standard but important part of reproducible computational biology documentation.