Bio.SeqIO 套件

子模組

模組內容

以 SeqRecord 物件進行序列輸入/輸出。

Bio.SeqIO 的文件也位於SeqIO,並且在我們的教學中有整章說明

輸入

主要函式為 Bio.SeqIO.parse(…),它接受一個輸入檔案控制代碼 (或在最近版本的 Biopython 中,也可以接受一個檔名作為字串) 和格式字串。這會返回一個迭代器,產生 SeqRecord 物件。

>>> from Bio import SeqIO
>>> for record in SeqIO.parse("Fasta/f002", "fasta"):
...     print("%s %i" % (record.id, len(record)))
gi|1348912|gb|G26680|G26680 633
gi|1348917|gb|G26685|G26685 413
gi|1592936|gb|G29385|G29385 471

請注意,parse() 函式會使用預設設定,為該格式調用相關的剖析器。您可能需要更多控制,在這種情況下,您需要直接建立一個特定格式的序列迭代器。

這些剖析器中有些是低階剖析器的包裝函式,用於建立 SeqRecord 物件,以實現一致的 SeqIO 介面。在執行時間至關重要的情況下 (例如大型 FASTA 或 FASTQ 檔案),調用這些底層剖析器會快得多 - 在這種情況下,這些產生器函式會返回字串元組

>>> from Bio.SeqIO.FastaIO import SimpleFastaParser
>>> from Bio.SeqIO.QualityIO import FastqGeneralIterator

輸入 - 單一記錄

如果您預期您的檔案只包含一個記錄,那麼我們會提供以下「輔助」函式,它會返回一個單一的 SeqRecord,如果沒有記錄或有多個記錄,則會引發例外

>>> from Bio import SeqIO
>>> record = SeqIO.read("Fasta/f001", "fasta")
>>> print("%s %i" % (record.id, len(record)))
gi|3318709|pdb|1A91| 79

當您預期只有一個記錄時 (並且會將多個記錄視為錯誤),此樣式很有用。例如,當處理細菌基因體或染色體的 GenBank 檔案時,通常只有一個記錄。或者,當從網際網路下載單一記錄時,請將其與控制代碼搭配使用。

但是,如果您只是想要包含多個記錄的檔案中的第一個記錄,請在迭代器上使用 next() 函式

>>> from Bio import SeqIO
>>> record = next(SeqIO.parse("Fasta/f002", "fasta"))
>>> print("%s %i" % (record.id, len(record)))
gi|1348912|gb|G26680|G26680 633

只要檔案至少包含一個記錄,上述程式碼就會正常運作。請注意,如果有多個記錄,則會自動忽略剩餘的記錄。

輸入 - 多個記錄

對於具有多個記錄的非交錯檔案 (例如 Fasta、GenBank、EMBL),使用序列迭代器可以節省大量記憶體 (RAM)。對於交錯檔案格式 (例如大多數多重比對檔案格式) 的好處較少。但是,迭代器只允許您逐個存取記錄。

如果您想要依數字隨機存取記錄,請將其轉換為清單

>>> from Bio import SeqIO
>>> records = list(SeqIO.parse("Fasta/f002", "fasta"))
>>> len(records)
3
>>> print(records[1].id)
gi|1348917|gb|G26685|G26685

如果您想要依索引鍵 (例如記錄 ID) 隨機存取記錄,請將迭代器轉換為字典

>>> from Bio import SeqIO
>>> record_dict = SeqIO.to_dict(SeqIO.parse("Fasta/f002", "fasta"))
>>> len(record_dict)
3
>>> print(len(record_dict["gi|1348917|gb|G26685|G26685"]))
413

但是,使用 list() 或 to_dict() 函式會將所有記錄一次載入記憶體,因此對於非常大的檔案而言是不可能的。相反,對於某些檔案格式,Bio.SeqIO 提供了一種索引方法,可提供對任何記錄的字典式存取。例如:

>>> from Bio import SeqIO
>>> record_dict = SeqIO.index("Fasta/f002", "fasta")
>>> len(record_dict)
3
>>> print(len(record_dict["gi|1348917|gb|G26685|G26685"]))
413
>>> record_dict.close()

許多(但非全部)支援的輸入檔案格式都可以像這樣建立索引。例如,「fasta」、「fastq」、「qual」甚至二進位格式「sff」都可以,但像「phylip」、「clustalw」和「nexus」這類比對格式則不行。

在大多數情況下,您也可以使用 SeqIO.index 從檔案中取得記錄作為原始字串(而非 SeqRecord)。這在以下情況會很有用:例如,從 SeqIO 無法輸出檔案格式的檔案(例如,純文字 SwissProt 格式「swiss」)中提取記錄的子集,或者在需要保持輸出與輸入 100% 相同的情況下。例如,

>>> from Bio import SeqIO
>>> record_dict = SeqIO.index("Fasta/f002", "fasta")
>>> len(record_dict)
3
>>> print(record_dict.get_raw("gi|1348917|gb|G26685|G26685").decode())
>gi|1348917|gb|G26685|G26685 human STS STS_D11734.
CGGAGCCAGCGAGCATATGCTGCATGAGGACCTTTCTATCTTACATTATGGCTGGGAATCTTACTCTTTC
ATCTGATACCTTGTTCAGATTTCAAAATAGTTGTAGCCTTATCCTGGTTTTACAGATGTGAAACTTTCAA
GAGATTTACTGACTTTCCTAGAATAGTTTCTCTACTGGAAACCTGATGCTTTTATAAGCCATTGTGATTA
GGATGACTGTTACAGGCTTAGCTTTGTGTGAAANCCAGTCACCTTTCTCCTAGGTAATGAGTAGTGCTGT
TCATATTACTNTAAGTTCTATAGCATACTTGCNATCCTTTANCCATGCTTATCATANGTACCATTTGAGG
AATTGNTTTGCCCTTTTGGGTTTNTTNTTGGTAAANNNTTCCCGGGTGGGGGNGGTNNNGAAA

>>> print(record_dict["gi|1348917|gb|G26685|G26685"].format("fasta"))
>gi|1348917|gb|G26685|G26685 human STS STS_D11734.
CGGAGCCAGCGAGCATATGCTGCATGAGGACCTTTCTATCTTACATTATGGCTGGGAATC
TTACTCTTTCATCTGATACCTTGTTCAGATTTCAAAATAGTTGTAGCCTTATCCTGGTTT
TACAGATGTGAAACTTTCAAGAGATTTACTGACTTTCCTAGAATAGTTTCTCTACTGGAA
ACCTGATGCTTTTATAAGCCATTGTGATTAGGATGACTGTTACAGGCTTAGCTTTGTGTG
AAANCCAGTCACCTTTCTCCTAGGTAATGAGTAGTGCTGTTCATATTACTNTAAGTTCTA
TAGCATACTTGCNATCCTTTANCCATGCTTATCATANGTACCATTTGAGGAATTGNTTTG
CCCTTTTGGGTTTNTTNTTGGTAAANNNTTCCCGGGTGGGGGNGGTNNNGAAA

>>> record_dict.close()

這裡,原始檔案和 Biopython 的輸出在換行上有所不同。另請注意,get_raw 方法將返回一個位元組物件,因此需要使用 decode 將其轉換為字串。

另請注意,get_raw 方法會保留換行符號。此範例 FASTQ 檔案使用 Unix 風格的結尾(僅 b"n"),

>>> from Bio import SeqIO
>>> fastq_dict = SeqIO.index("Quality/example.fastq", "fastq")
>>> len(fastq_dict)
3
>>> raw = fastq_dict.get_raw("EAS54_6_R1_2_1_540_792")
>>> raw.count(b"\n")
4
>>> raw.count(b"\r\n")
0
>>> b"\r" in raw
False
>>> len(raw)
78
>>> fastq_dict.close()

這是同一個檔案,但使用 DOS/Windows 的換行符號(b"rn"),

>>> from Bio import SeqIO
>>> fastq_dict = SeqIO.index("Quality/example_dos.fastq", "fastq")
>>> len(fastq_dict)
3
>>> raw = fastq_dict.get_raw("EAS54_6_R1_2_1_540_792")
>>> raw.count(b"\n")
4
>>> raw.count(b"\r\n")
4
>>> b"\r\n" in raw
True
>>> len(raw)
82
>>> fastq_dict.close()

由於每個換行符號使用兩個位元組,因此檔案會比只有一個位元組的 Unix 對等檔案更長。

輸入 - 比對

您可以使用 Bio.AlignIO 將比對檔案讀取為比對物件。或者,透過 Bio.SeqIO 讀取比對檔案格式會為每個比對的每一列提供一個 SeqRecord。

>>> from Bio import SeqIO
>>> for record in SeqIO.parse("Clustalw/hedgehog.aln", "clustal"):
...     print("%s %i" % (record.id, len(record)))
gi|167877390|gb|EDS40773.1| 447
gi|167234445|ref|NP_001107837. 447
gi|74100009|gb|AAZ99217.1| 447
gi|13990994|dbj|BAA33523.2| 447
gi|56122354|gb|AAV74328.1| 447

輸出

使用函數 Bio.SeqIO.write(…),它會接收一組完整的 SeqRecord 物件(以列表或迭代器的形式)、一個輸出檔案控制代碼(或在較新版本的 Biopython 中,以字串形式表示的輸出檔案名稱)以及檔案格式

from Bio import SeqIO
records = ...
SeqIO.write(records, "example.faa", "fasta")

或者,使用控制代碼

from Bio import SeqIO
records = ...
with open("example.faa", "w") as handle:
  SeqIO.write(records, handle, "fasta")

您應該呼叫此函數一次(帶上所有記錄),如果使用控制代碼,請務必關閉它以將資料刷新到硬碟。

輸出 - 進階

在單一檔案上多次呼叫 write() 的效果會因檔案格式而異,最好避免這樣做,除非您有充分的理由。

如果您提供檔案名稱,則每次呼叫 write() 時都會覆寫現有檔案。對於循序檔案格式(例如 fasta、genbank),每個「記錄區塊」都保存單一序列。對於這些檔案,重複使用相同的控制代碼多次呼叫 write() 可能會是安全的。

但是,嘗試對某些比對格式(例如 phylip、clustal、stockholm)執行此操作,會導致多個多序列比對串連在一起。此類檔案是由 PHYLIP 程式套件為引導分析而建立的,但透過 Bio.AlignIO 執行此操作會更清楚。

更糟糕的是,許多檔案格式都有明確的標頭和/或頁尾結構(例如,任何 XMl 格式,以及大多數二進位檔案格式,如 SFF)。在此處多次呼叫 write() 將導致檔案無效。

轉換

Bio.SeqIO.convert(…) 函數為簡單的檔案格式轉換提供了簡便的介面。此外,它還可能使用特定於檔案格式的優化,因此這也應該是最快的方法。

但一般而言,您可以將 Bio.SeqIO.parse(…) 函數與 Bio.SeqIO.write(…) 函數結合使用,進行序列檔案轉換。使用產生器表達式或產生器函數可提供一種記憶體效率高的方法,以在流程中執行篩選或其他額外操作。

檔案格式

在指定檔案格式時,請使用小寫字串。相同的格式名稱也用於 Bio.AlignIO 中,包括以下格式

  • abi - Applied Biosystem 的定序追蹤格式

  • abi-trim - 與「abi」相同,但使用 Mott 演算法進行品質修剪

  • ace - 從 ACE 組裝檔案讀取重疊群序列。

  • cif-atom - 使用 Bio.PDB.MMCIFParser 來根據原子座標確定結構中出現的(部分)蛋白質序列。

  • cif-seqres - 讀取巨分子結晶學資訊檔案 (mmCIF) 以確定 _pdbx_poly_seq_scheme 記錄定義的完整蛋白質序列。

  • embl - EMBL 平面檔案格式。在內部使用 Bio.GenBank。

  • fasta - 一般序列檔案格式,其中每條記錄都以「>」字元開頭的識別碼行開始,後跟序列行。

  • fasta-2line - 使用每條記錄恰好兩行的 FASTA 格式的更嚴格解釋(不換行)。

  • fastq - Sanger 使用的「類似 FASTA」格式,也會儲存 PHRED 序列品質值(ASCII 偏移量為 33)。

  • fastq-sanger - 與 BioPerl 和 EMBOSS 一致的「fastq」別名

  • fastq-solexa - FASTQ 格式的原始 Solexa/Illumnia 變體,它使用 ASCII 偏移量 64 編碼 Solexa 品質分數(不是 PHRED 品質分數)。

  • fastq-illumina - FASTQ 格式的 Solexa/Illumina 1.3 至 1.7 變體,它使用 ASCII 偏移量 64(不是 33)編碼 PHRED 品質分數。請注意,自 CASAVA 管線的 1.8 版起,Illumina 將使用標準 Sanger 編碼產生 FASTQ 檔案。

  • gck - Gene Construction Kit 的格式。

  • genbank - GenBank 或 GenPept 平面檔案格式。

  • gb - 與 NCBI Entrez Utilities 一致的「genbank」別名

  • gfa1 - 圖形片段組裝 1.x 版。僅剖析片段行,並忽略所有連結資訊。

  • gfa2 - 圖形片段組裝 2.0 版。僅剖析片段行,並忽略所有連結資訊。

  • ig - IntelliGenetics 檔案格式,顯然與 MASE 比對格式相同。

  • imgt - 來自 IMGT 的 EMBL 類似格式,其中特徵表縮排更多,以允許更長的功能類型。

  • nib - UCSC 的核苷酸序列 nib 檔案格式,它使用一個半位元組(4 位元)來表示每個核苷酸,並在一個位元組中儲存兩個核苷酸。

  • pdb-seqres - 讀取蛋白質資料庫 (PDB) 檔案以確定標頭中顯示的完整蛋白質序列(無相依性)。

  • pdb-atom - 使用 Bio.PDB 根據檔案的原子座標區段中顯示的結構來確定(部分)蛋白質序列(Bio.PDB 需要 NumPy)。

  • phd - PHRED 的輸出,由 PHRAP 和 CONSED 用於輸入。

  • pir - 國家生物醫學研究基金會 (NBRF) 針對蛋白質資訊資源 (PIR) 資料庫(現在是 UniProt 的一部分)引入的「類似 FASTA」格式。

  • seqxml - Schmitt 等人 (2011) 中描述的 SeqXML 簡單 XML 格式。

  • sff - 標準流程圖格式 (SFF),Roche 454 的典型輸出。

  • sff-trim - 套用給定修剪的標準流程圖格式 (SFF)。

  • snapgene - SnapGene 的原生格式。

  • swiss - 純文字 Swiss-Prot 又名 UniProt 格式。

  • tab - 簡單的兩欄以 tab 分隔的序列檔案,其中每一行都保存記錄的識別碼和序列。例如,當以最小的 tab 分隔文字檔案儲存微陣列探針時,Aligent 的 eArray 軟體會使用此檔案。

  • qual - 一種「類似 FASTA」的格式,保存來自定序 DNA 的 PHRED 品質值,但沒有實際序列(通常在單獨的 FASTA 檔案中提供)。

  • uniprot-xml - UniProt XML 格式(取代我們稱為「swiss」的 SwissProt 純文字格式)

  • xdna - DNA Strider 和 SerialCloner 的原生格式。

請注意,雖然 Bio.SeqIO 可以讀取上述所有檔案格式,但它無法寫入所有格式。

您也可以使用 Bio.AlignIO 支援的任何檔案格式,例如「nexus」、「phylip」和「stockholm」,這可讓您存取構成每個比對的個別序列作為 SeqRecord。

Bio.SeqIO.write(sequences: Iterable[SeqRecord] | SeqRecord, handle: IO[str] | PathLike | str | bytes, format: str) int

將完整的序列集寫入檔案。

引數
  • sequences - SeqRecord 物件的列表(或迭代器),或單一 SeqRecord。

  • handle - 要寫入的檔案控制代碼物件,或作為字串的檔案名稱。

  • format - 描述要寫入的檔案格式的小寫字串。

請注意,如果提供檔案控制代碼,您的程式碼應在呼叫此函數後關閉控制代碼(以確保資料刷新到磁碟)。

傳回寫入的記錄數(為整數)。

Bio.SeqIO.parse(handle, format, alphabet=None)

將序列檔案轉換為傳回 SeqRecord 的迭代器。

引數
  • handle - 檔案的控制代碼,或作為字串的檔案名稱(請注意,舊版本的 Biopython 僅採用控制代碼)。

  • format - 描述檔案格式的小寫字串。

  • alphabet - 不再使用,應為 None。

典型用法,開啟檔案以讀取,並迴圈處理記錄

>>> from Bio import SeqIO
>>> filename = "Fasta/sweetpea.nu"
>>> for record in SeqIO.parse(filename, "fasta"):
...    print("ID %s" % record.id)
...    print("Sequence length %i" % len(record))
ID gi|3176602|gb|U78617.1|LOU78617
Sequence length 309

對於如 twobit 這類僅在需要時才讀取檔案內容的延遲載入檔案格式,請確保在擷取序列資料時檔案保持開啟。

如果您有一個包含檔案內容的字串 ‘data’,您必須先將其轉換為一個 handle 才能進行解析。

>>> data = ">Alpha\nACCGGATGTA\n>Beta\nAGGCTCGGTTA\n"
>>> from Bio import SeqIO
>>> from io import StringIO
>>> for record in SeqIO.parse(StringIO(data), "fasta"):
...     print("%s %s" % (record.id, record.seq))
Alpha ACCGGATGTA
Beta AGGCTCGGTTA

當您預期只有單一記錄時,請使用 Bio.SeqIO.read(…) 函數。

Bio.SeqIO.read(handle, format, alphabet=None)

將序列檔案轉換為單一 SeqRecord。

引數
  • handle - 檔案的控制代碼,或作為字串的檔案名稱(請注意,舊版本的 Biopython 僅採用控制代碼)。

  • format - 描述檔案格式的字串。

  • alphabet - 不再使用,應為 None。

此函數用於解析包含正好一筆記錄的序列檔案。例如,讀取 GenBank 檔案。

>>> from Bio import SeqIO
>>> record = SeqIO.read("GenBank/arab1.gb", "genbank")
>>> print("ID %s" % record.id)
ID AC007323.5
>>> print("Sequence length %i" % len(record))
Sequence length 86436

如果 handle 不包含任何記錄,或包含超過一筆記錄,則會引發例外。例如:

>>> from Bio import SeqIO
>>> record = SeqIO.read("GenBank/cor6_6.gb", "genbank")
Traceback (most recent call last):
    ...
ValueError: More than one record found in handle

然而,如果您想要從包含多筆記錄的檔案中取得第一筆記錄,此函數會引發例外(如上例所示)。請改用:

>>> from Bio import SeqIO
>>> record = next(SeqIO.parse("GenBank/cor6_6.gb", "genbank"))
>>> print("First record's ID %s" % record.id)
First record's ID X55053.1

如果您想從 handle 讀取多筆記錄,請使用 Bio.SeqIO.parse(handle, format) 函數。

Bio.SeqIO.to_dict(sequences, key_function=None)

將序列迭代器或列表轉換為字典。

引數
  • sequences - 一個返回 SeqRecord 物件的迭代器,或只是一個 SeqRecord 物件的列表。

  • key_function - 可選的回調函數,當給定 SeqRecord 時,應返回字典的唯一鍵。

例如:key_function = lambda rec : rec.name 或 key_function = lambda rec : rec.description.split()[0]

如果省略 key_function,則會使用 record.id,假設返回的記錄物件是具有唯一 id 的 SeqRecords。

如果有重複的鍵,則會引發錯誤。

自 Python 3.7 起,預設的 dict 類別會維護鍵的順序,表示此字典將反映提供給它的記錄順序。對於 CPython 和 PyPy,這已在 Python 3.6 中實作,因此實際上您可以始終假設記錄順序會被保留。

範例用法,預設使用 record.id 作為鍵:

>>> from Bio import SeqIO
>>> filename = "GenBank/cor6_6.gb"
>>> format = "genbank"
>>> id_dict = SeqIO.to_dict(SeqIO.parse(filename, format))
>>> print(list(id_dict))
['X55053.1', 'X62281.1', 'M81224.1', 'AJ237582.1', 'L31939.1', 'AF297471.1']
>>> print(id_dict["L31939.1"].description)
Brassica rapa (clone bif72) kin mRNA, complete cds

一個更複雜的例子,使用 key_function 參數以使用序列校驗和作為字典的鍵:

>>> from Bio import SeqIO
>>> from Bio.SeqUtils.CheckSum import seguid
>>> filename = "GenBank/cor6_6.gb"
>>> format = "genbank"
>>> seguid_dict = SeqIO.to_dict(SeqIO.parse(filename, format),
...               key_function = lambda rec : seguid(rec.seq))
>>> for key, record in sorted(seguid_dict.items()):
...     print("%s %s" % (key, record.id))
/wQvmrl87QWcm9llO4/efg23Vgg AJ237582.1
BUg6YxXSKWEcFFH0L08JzaLGhQs L31939.1
SabZaA4V2eLE9/2Fm5FnyYy07J4 X55053.1
TtWsXo45S3ZclIBy4X/WJc39+CY M81224.1
l7gjJFE6W/S1jJn5+1ASrUKW/FA X62281.1
uVEYeAQSV5EDQOnFoeMmVea+Oow AF297471.1

此方法不適用於非常大的序列集合,因為所有 SeqRecord 物件都保存在記憶體中。請考慮改用 Bio.SeqIO.index() 函數(如果它支援您的特定檔案格式)。

此字典將反映提供給它的記錄順序。

Bio.SeqIO.index(filename, format, alphabet=None, key_function=None)

索引一個序列檔案並返回一個類似字典的物件。

引數
  • filename - 提供要索引的檔案名稱的字串

  • format - 描述檔案格式的小寫字串

  • alphabet - 不再使用,請保留為 None

  • key_function - 可選的回調函數,當給定 SeqRecord 識別符字串時,應返回字典的唯一鍵。

此索引函數將返回一個類似字典的物件,將 SeqRecord 物件作為值。

從 Biopython 1.69 開始,當迭代條目時,這將保留檔案中記錄的順序。

>>> from Bio import SeqIO
>>> records = SeqIO.index("Quality/example.fastq", "fastq")
>>> len(records)
3
>>> list(records)  # make a list of the keys
['EAS54_6_R1_2_1_413_324', 'EAS54_6_R1_2_1_540_792', 'EAS54_6_R1_2_1_443_348']
>>> print(records["EAS54_6_R1_2_1_540_792"].format("fasta"))
>EAS54_6_R1_2_1_540_792
TTGGCAGGCCAAGGCCGATGGATCA

>>> "EAS54_6_R1_2_1_540_792" in records
True
>>> print(records.get("Missing", None))
None
>>> records.close()

如果檔案是 BGZF 壓縮的,則會自動偵測到。不支援普通的 GZIP 檔案。

>>> from Bio import SeqIO
>>> records = SeqIO.index("Quality/example.fastq.bgz", "fastq")
>>> len(records)
3
>>> print(records["EAS54_6_R1_2_1_540_792"].seq)
TTGGCAGGCCAAGGCCGATGGATCA
>>> records.close()

當您呼叫索引函數時,它會掃描整個檔案,記錄每個記錄的位置。當您透過字典方法存取特定記錄時,程式碼會跳到檔案的適當部分,然後將該部分解析為 SeqRecord。

請注意,並非 Bio.SeqIO 支援的所有輸入格式都可與此索引函數一起使用。它設計為僅適用於循序檔案格式(例如,“fasta”、“gb”、“fastq”),而不適用於任何交錯檔案格式(例如,對齊格式,如“clustal”)。

對於小檔案,使用記憶體中的 Python 字典可能更有效率,例如:

>>> from Bio import SeqIO
>>> records = SeqIO.to_dict(SeqIO.parse("Quality/example.fastq", "fastq"))
>>> len(records)
3
>>> list(records)  # make a list of the keys
['EAS54_6_R1_2_1_413_324', 'EAS54_6_R1_2_1_540_792', 'EAS54_6_R1_2_1_443_348']
>>> print(records["EAS54_6_R1_2_1_540_792"].format("fasta"))
>EAS54_6_R1_2_1_540_792
TTGGCAGGCCAAGGCCGATGGATCA

與 to_dict() 函數一樣,預設情況下,每個記錄的 id 字串會用作鍵。您可以指定一個回調函數來將此(記錄識別符字串)轉換為您偏好的鍵。例如:

>>> from Bio import SeqIO
>>> def make_tuple(identifier):
...     parts = identifier.split("_")
...     return int(parts[-2]), int(parts[-1])
>>> records = SeqIO.index("Quality/example.fastq", "fastq",
...                       key_function=make_tuple)
>>> len(records)
3
>>> list(records)  # make a list of the keys
[(413, 324), (540, 792), (443, 348)]
>>> print(records[(540, 792)].format("fasta"))
>EAS54_6_R1_2_1_540_792
TTGGCAGGCCAAGGCCGATGGATCA

>>> (540, 792) in records
True
>>> "EAS54_6_R1_2_1_540_792" in records
False
>>> print(records.get("Missing", None))
None
>>> records.close()

另一個常見的用例是索引 NCBI 樣式的 FASTA 檔案,您可能希望從 FASTA 識別符中擷取 GI 編號以用作字典的鍵。

請注意,與 to_dict() 函數不同,這裡的 key_function 不會獲得完整的 SeqRecord 來用於產生鍵。這樣做會造成嚴重的效能損失,因為它需要在建立索引時完全解析檔案。目前通常會避免這種情況。

另請參閱:Bio.SeqIO.index_db() 和 Bio.SeqIO.to_dict()

Bio.SeqIO.index_db(index_filename, filenames=None, format=None, alphabet=None, key_function=None)

索引多個序列檔案並返回一個類似字典的物件。

索引儲存在 SQLite 資料庫中,而不是在記憶體中(如 Bio.SeqIO.index(…) 函數中)。

引數
  • index_filename - 儲存 SQLite 索引的位置

  • filenames - 指定要索引的檔案清單的字串,或者在索引單一檔案時,可以將其作為字串提供。(如果重新載入現有的索引則為選用,但必須匹配)

  • format - 描述檔案格式的小寫字串(如果重新載入現有的索引則為選用,但必須匹配)

  • alphabet - 不再使用,請保留為 None。

  • key_function - 可選的回調函數,當給定 SeqRecord 識別符字串時,應返回字典的唯一鍵。

此索引函數將返回一個類似字典的物件,將 SeqRecord 物件作為值

>>> from Bio import SeqIO
>>> files = ["GenBank/NC_000932.faa", "GenBank/NC_005816.faa"]
>>> def get_gi(name):
...     parts = name.split("|")
...     i = parts.index("gi")
...     assert i != -1
...     return parts[i+1]
>>> idx_name = ":memory:" #use an in memory SQLite DB for this test
>>> records = SeqIO.index_db(idx_name, files, "fasta", key_function=get_gi)
>>> len(records)
95
>>> records["7525076"].description
'gi|7525076|ref|NP_051101.1| Ycf2 [Arabidopsis thaliana]'
>>> records["45478717"].description
'gi|45478717|ref|NP_995572.1| pesticin [Yersinia pestis biovar Microtus str. 91001]'
>>> records.close()

在此範例中,兩個檔案分別包含 85 和 10 筆記錄。

支援 BGZF 壓縮檔案,並且會自動偵測。不支援普通的 GZIP 壓縮檔案。

另請參閱:Bio.SeqIO.index() 和 Bio.SeqIO.to_dict(),以及 Python 模組 glob,該模組對於建立檔案清單很有用。

Bio.SeqIO.convert(in_file, in_format, out_file, out_format, molecule_type=None)

在兩種序列檔案格式之間轉換,返回記錄數。

引數
  • in_file - 輸入 handle 或檔案名稱

  • in_format - 輸入檔案格式,小寫字串

  • out_file - 輸出 handle 或檔案名稱

  • out_format - 輸出檔案格式,小寫字串

  • molecule_type - 要套用的可選分子類型,包含 “DNA”、“RNA” 或 “protein” 的字串。

注意 - 如果您提供輸出檔案名稱,它將被開啟,這將覆蓋任何現有的檔案,恕不另行警告。

這裡的想法是,雖然這樣做會有效

from Bio import SeqIO
records = SeqIO.parse(in_handle, in_format)
count = SeqIO.write(records, out_handle, out_format)

但寫成這樣會更簡短

from Bio import SeqIO
count = SeqIO.convert(in_handle, in_format, out_handle, out_format)

此外,Bio.SeqIO.convert 對於某些轉換速度更快,因為它可以進行一些最佳化。

例如,從檔案名稱轉換為 handle

>>> from Bio import SeqIO
>>> from io import StringIO
>>> handle = StringIO("")
>>> SeqIO.convert("Quality/example.fastq", "fastq", handle, "fasta")
3
>>> print(handle.getvalue())
>EAS54_6_R1_2_1_413_324
CCCTTCTTGTCTTCAGCGTTTCTCC
>EAS54_6_R1_2_1_540_792
TTGGCAGGCCAAGGCCGATGGATCA
>EAS54_6_R1_2_1_443_348
GTTGCTTCTGGCGTGGGTGGGGGGG

請注意,某些格式(如 SeqXML)要求您在剖析器無法確定時指定分子類型

>>> from Bio import SeqIO
>>> from io import BytesIO
>>> handle = BytesIO()
>>> SeqIO.convert("Quality/example.fastq", "fastq", handle, "seqxml", "DNA")
3