Bio.PDB.internal_coords 模組

支援蛋白質結構內部坐標的類別。

內部坐標包括沿著蛋白質骨架的 Psi、Omega 和 Phi 二面角,沿著側鏈的 Chi 角,以及定義蛋白質鏈的所有 3 原子角和鍵長。這些例程可以從原子 XYZ 坐標計算內部坐標,並從內部坐標計算原子 XYZ 坐標。

次要優點包括能夠比對和比較 3D 結構中的殘基環境、支援 2D 原子距離圖、將距離圖加上手性資訊轉換為結構、產生用於 3D 列印的結構 OpenSCAD 描述,以及以內部坐標資料檔案讀取/寫入結構。

用法

from Bio.PDB.PDBParser import PDBParser
from Bio.PDB.Chain import Chain
from Bio.PDB.internal_coords import *
from Bio.PDB.PICIO import write_PIC, read_PIC, read_PIC_seq
from Bio.PDB.ic_rebuild import write_PDB, IC_duplicate, structure_rebuild_test
from Bio.PDB.SCADIO import write_SCAD
from Bio.Seq import Seq
from Bio.SeqRecord import SeqRecord
from Bio.PDB.PDBIO import PDBIO
import numpy as np

# load a structure as normal, get first chain
parser = PDBParser()
myProtein = parser.get_structure("7rsa", "pdb7rsa.ent")
myChain = myProtein[0]["A"]

# compute bond lengths, angles, dihedral angles
myChain.atom_to_internal_coordinates(verbose=True)

# check myChain makes sense (can get angles and rebuild same structure)
resultDict = structure_rebuild_test(myChain)
assert resultDict['pass'] == True

# get residue 1 chi2 angle
r1 = next(myChain.get_residues())
r1chi2 = r1.internal_coord.get_angle("chi2")

# rotate residue 1 chi2 angle by 120 degrees (loops w/in +/-180)
r1.internal_coord.set_angle("chi2", r1chi2 + 120.0)
# or
r1.internal_coord.bond_rotate("chi2", 120.0)
# update myChain XYZ coordinates with chi2 changed
myChain.internal_to_atom_coordinates()
# write new conformation with PDBIO
write_PDB(myProtein, "myChain.pdb")
# or just the ATOM records without headers:
io = PDBIO()
io.set_structure(myProtein)
io.save("myChain2.pdb")

# write chain as 'protein internal coordinates' (.pic) file
write_PIC(myProtein, "myChain.pic")
# read .pic file
myProtein2 = read_PIC("myChain.pic")

# create default structure for random sequence by reading as .pic file
myProtein3 = read_PIC_seq(
    SeqRecord(
        Seq("GAVLIMFPSTCNQYWDEHKR"),
        id="1RND",
        description="my random sequence",
    )
)
myProtein3.internal_to_atom_coordinates()
write_PDB(myProtein3, "myRandom.pdb")

# access the all-dihedrals array for the chain, e.g. residue 1 chi2 angle:
r1chi2_obj = r1.internal_coord.pick_angle("chi2")
# or same thing: r1chi2_obj = r1.internal_coord.pick_angle("CA:CB:CG:CD")
r1chi2_key = r1chi2_obj.atomkeys
# r1chi2_key is tuple of AtomKeys (1_K_CA, 1_K_CB, 1_K_CG, 1_K_CD)
r1chi2_index = myChain.internal_coord.dihedraNdx[r1chi2_key]
# or same thing: r1chi2_index = r1chi2_obj.ndx
r1chi2_value = myChain.internal_coord.dihedraAngle[r1chi2_index]
# also true: r1chi2_obj == myChain.internal_coord.dihedra[r1chi2_index]

# access the array of all atoms for the chain, e.g. residue 1 C-beta
r1_cBeta_index = myChain.internal_coord.atomArrayIndex[AtomKey("1_K_CB")]
r1_cBeta_coords = myChain.internal_coord.atomArray[r1_cBeta_index]
# r1_cBeta_coords = [ x, y, z, 1.0 ]

# the Biopython Atom coord array is now a view into atomArray, so
assert r1_cBeta_coords[1] == r1["CB"].coord[1]
r1_cBeta_coords[1] += 1.0  # change the Y coord 1 angstrom
assert r1_cBeta_coords[1] == r1["CB"].coord[1]
# they are always the same (they share the same memory)
r1_cBeta_coords[1] -= 1.0  # restore

# create a selector to filter just the C-alpha atoms from the all atom array
atmNameNdx = AtomKey.fields.atm
atomArrayIndex = myChain.internal_coord.atomArrayIndex
CaSelect = [
    atomArrayIndex.get(k) for k in atomArrayIndex.keys() if k.akl[atmNameNdx] == "CA"
]
# now the ordered array of C-alpha atom coordinates is:
CA_coords = myChain.internal_coord.atomArray[CaSelect]
# note this uses Numpy fancy indexing, so CA_coords is a new copy

# create a C-alpha distance plot
caDistances = myChain.internal_coord.distance_plot(CaSelect)
# display with e.g. MatPlotLib:
# import matplotlib.pyplot as plt
# plt.imshow(caDistances, cmap="hot", interpolation="nearest")
# plt.show()

# build structure from distance plot:
## create the all-atom distance plot
distances = myChain.internal_coord.distance_plot()
## get the sign of the dihedral angles
chirality = myChain.internal_coord.dihedral_signs()
## get new, empty data structure : copy data structure from myChain
myChain2 = IC_duplicate(myChain)[0]["A"]
cic2 = myChain2.internal_coord
## clear the new atomArray and di/hedra value arrays, just for proof
cic2.atomArray = np.zeros((cic2.AAsiz, 4), dtype=np.float64)
cic2.dihedraAngle[:] = 0.0
cic2.hedraAngle[:] = 0.0
cic2.hedraL12[:] = 0.0
cic2.hedraL23[:] = 0.0
## copy just the first N-Ca-C coords so structures will superimpose:
cic2.copy_initNCaCs(myChain.internal_coord)
## copy distances to chain arrays:
cic2.distplot_to_dh_arrays(distances, chirality)
## compute angles and dihedral angles from distances:
cic2.distance_to_internal_coordinates()
## generate XYZ coordinates from internal coordinates:
myChain2.internal_to_atom_coordinates()
## confirm result atomArray matches original structure:
assert np.allclose(cic2.atomArray, myChain.internal_coord.atomArray)

# superimpose all phe-phe pairs - quick hack just to demonstrate concept
# for analyzing pairwise residue interactions.  Generates PDB ATOM records
# placing each PHE at origin and showing all other PHEs in environment
## shorthand for key variables:
cic = myChain.internal_coord
resNameNdx = AtomKey.fields.resname
aaNdx = cic.atomArrayIndex
## select just PHE atoms:
pheAtomSelect = [aaNdx.get(k) for k in aaNdx.keys() if k.akl[resNameNdx] == "F"]
aaF = cic.atomArray[ pheAtomSelect ]  # numpy fancy indexing makes COPY not view

for ric in cic.ordered_aa_ic_list:  # internal_coords version of get_residues()
    if ric.rbase[2] == "F":  # if PHE, get transform matrices for chi1 dihedral
        chi1 = ric.pick_angle("N:CA:CB:CG")  # chi1 space has C-alpha at origin
        cst = np.transpose(chi1.cst)  # transform TO chi1 space
        # rcst = np.transpose(chi1.rcst)  # transform FROM chi1 space
        cic.atomArray[pheAtomSelect] = aaF.dot(cst)  # transform just the PHEs
        for res in myChain.get_residues():  # print PHEs in new coordinate space
            if res.resname in ["PHE"]:
                print(res.internal_coord.pdb_residue_string())
        cic.atomArray[pheAtomSelect] = aaF  # restore coordinate space from copy

# write OpenSCAD program of spheres and cylinders to 3d print myChain backbone
## set atom load filter to accept backbone only:
IC_Residue.accept_atoms = IC_Residue.accept_backbone
## delete existing data to force re-read of all atoms:
myChain.internal_coord = None
write_SCAD(myChain, "myChain.scad", scale=10.0)

請參閱《Biopython 教學與手冊》的「內部坐標模組」部分,以了解進一步的討論。

術語和關鍵資料結構: 內部坐標定義在跨越殘基或沿著側鏈遵循公認命名法的原子序列上。為了管理這些序列並支援 Biopython 的無序機制,實作了 AtomKey 規範器,以在單一物件中捕獲殘基、原子和變體識別資訊。一個 Hedron 物件被指定為三個連續的 AtomKey,包含兩個鍵長以及它們之間的鍵角。一個 Dihedron 由四個連續的 AtomKey 組成,將兩個 Hedra 與它們之間的二面角連結起來。

演算法概述: 內部坐標模組結合了 ic_data 檔案中連接原子的規範,作為 hedra 和 dihedra,並在此處提供例程,以在局部坐標系統與例如 PDB 或 mmCif 資料檔案中提供的全局坐標之間轉換這些原子集的 XYZ 坐標。局部坐標系統將 hedron 的中心原子放置在原點 (0,0,0),一個腳位於 +Z 軸上,另一個腳位於 XZ 平面上(請參閱 Hedron)。在局部坐標空間中測量和建立或操作 hedra 和 dihedra 很簡單,並且計算出的轉換矩陣可以從初始 N-Ca-C 原子的供應(PDB)坐標開始將這些子單元組裝成蛋白質鏈。

Psi 和 Phi 角定義在蛋白質鏈中相鄰殘基的原子上,例如,請參閱 pick_angle()ic_data,以了解殘基和骨架二面角之間的相關映射。

可以透過 IC_Chain.dCoordSpaceDihedron 屬性 .cst 和 .rcst 存取到和從上述二面體局部坐標空間的轉換,並且可以在殘基及其環境的比對和比較中應用,其程式碼如下:

chi1 = ric0.pick_angle("chi1") # chi1 space defined with CA at origin
cst = np.transpose(chi1.cst) # transform TO chi1 local space
newAtomCoords = oldAtomCoords.dot(cst)

核心演算法是在 1993-4 年間獨立開發的,用於 《蛋白質中胺基酸環境的三維描述的開發和應用》,Miller、Douthart 和 Dunker,《分子生物資訊學進展》,IOS Press,1994 年,ISBN 90 5199 172 x,第 9-30 頁。

定義了蛋白質內部坐標 (.pic) 檔案格式,以捕獲足夠的細節,從鏈起點坐標(第一個殘基 N、Ca、C XYZ 坐標)和剩餘的內部坐標重現 PDB 檔案。這些檔案在內部用於驗證是否可以從其內部坐標重新產生給定的結構。請參閱 PICIO,以了解如何讀取和寫入 .pic 檔案,並參閱 structure_rebuild_test(),以判斷特定的 PDB 或 mmCif 資料檔案是否具有足夠的資訊來在笛卡爾坐標和內部坐標之間相互轉換。

內部坐標也可以匯出為 OpenSCAD 資料陣列,用於產生 3D 列印蛋白質模型。提供 OpenSCAD 軟體作為產生此類模型的起點和概念驗證。請參閱 SCADIO 和這個 Thingiverse 專案,以取得更進階的範例。

請參閱 distance_plot()distance_to_internal_coordinates(),以了解如何將結構資料轉換為 2D 距離圖或從 2D 距離圖轉換結構資料。

下列類別包含用於處理內部坐標的核心功能,並且彼此之間具有足夠的關聯性和耦合性,因此將它們放在此模組中

IC_Chain:在 .internal_coord 屬性上擴充 Biopython 鏈。

管理連接的殘基序列和鏈斷裂;保存所有原子坐標和鍵幾何結構的 numpy 陣列。為了進行「平行」處理,IC_Chain 方法使用單一 numpy 命令來操作這些陣列。

IC_Residue:在 .internal_coord 屬性上擴充 Biopython 殘基。

存取內部坐標的每個殘基視圖,以及用於串行(殘基接殘基)組裝的方法。

Dihedron:形成二面角的四個連接原子。

二面角、局部坐標空間中的均勻原子坐標、對相關 Hedra 和 IC_Residue 的參照。殘基二面角、鍵角和鍵長的 Getter 方法。

Hedron:形成平面的三個連接原子。

包含局部坐標空間中的均勻原子坐標,以及鍵長和它們之間的夾角。

Edron:Hedron 和 Dihedron 類別的基底類別。

AtomKey 的元組,包含子項、字串 ID、主鏈成員資格布林值,以及 Hedra 和 Dihedra 共用的其他例程。實作豐富的比較。

AtomKey:用於參照原子序列的鍵(字典和字串)。

捕獲殘基和無序/佔用資訊,為 .pic 檔案提供無空白鍵,並實作豐富的比較。

自訂例外類別:HedronMatchErrorMissingAtomError

class Bio.PDB.internal_coords.IC_Chain(parent, verbose: bool = False)

基底類別:object

使用內部坐標資料擴充 Biopython 鏈的類別。

屬性:
chain: 物件參照

此物件繼承自 Biopython 的 Bio.PDB.Chain.Chain 物件

MaxPeptideBond: float

用於偵測鏈斷裂的類別屬性。若鏈為完全連續,但具有非常長的鍵(例如,用於 3D 列印 (OpenSCAD 輸出) 缺少殘基的結構),則可覆寫此屬性。MaxPeptideBond

ParallelAssembleResidues: bool

影響 internal_to_atom_coords 的類別屬性。較短的鏈(50 個殘基或更少)在組裝時,若不建立 numpy 陣列,速度會更快。且一次處理一個殘基的演算法也更容易理解和追蹤。清除此標誌(設定為 False)將會切換為序列演算法

ordered_aa_ic_list: list

內部座標演算法可以處理的 IC_Residue 物件(例如,不包含水)

initNCaC: N、Ca、C 原子鍵 (AtomKey) 元組 (NCaCKeys) 的列表。

NCaCKeys 會開始鏈片段(第一個殘基或在鏈斷裂之後)。這 3 個原子定義了連續鏈片段的座標空間,此空間最初由 PDB 或 mmCIF 檔案指定。

AAsiz = int

原子陣列的大小,此鏈中的原子數量

atomArray: numpy array

鏈中每個原子的均勻原子座標 ([x,, y, z, 1.0])

atomArrayIndex: dict

將 AtomKeys 對應至 atomArray 索引

hedra: dict

此鏈中形成多面體 (Hedra) 的殘基;以 AtomKeys 的 3 元組作為索引。

hedraLen: int

hedra 字典的長度

hedraNdx: dict

將 hedra AtomKeys 對應至 hedra 資料陣列中的數值索引,例如下方的 hedraL12

a2ha_map: [hedraLen x 3]

以 hedraNdx 順序排列的原子索引

dihedra: dict

此鏈中形成二面體 (Dihedra) 的殘基;以 AtomKeys 的 4 元組作為索引。

dihedraLen: int

dihedra 字典的長度

dihedraNdx: dict

將 dihedra AtomKeys 對應至 dihedra 資料陣列,例如 dihedraAngle

a2da_map[dihedraLen x 4]

以 dihedraNdx 順序排列的 AtomNdx

d2a_map[dihedraLen x [4]]

每個二面體的 AtomNdx(重新塑形的 a2da_map)

用於向量處理鏈二面體/多面體的 Numpy 陣列
hedraL12: numpy array

多面體第 1 個和第 2 個原子之間的鍵長

hedraAngle: numpy array

每個多面體的鍵角,以度為單位

hedraL23: numpy array

多面體第 2 個和第 3 個原子之間的鍵長

id3_dh_index: dict

將多面體鍵對應至以多面體開始的二面體列表,由 assemble 和 bond_rotate 用來尋找具有 h1 鍵的二面體

id32_dh_index: dict

類似 id3_dh_index,從 h2 鍵尋找二面體

hAtoms: numpy array

多面體的均勻原子座標 (3x4),中心原子位於原點

hAtomsR: numpy array

反向排列的 hAtoms

hAtoms_needs_update: numpy array of bool

指出 hAtoms 是否代表 hedraL12/A/L23

dihedraAngle: numpy array

每個二面體的二面角(度)

dAtoms: numpy array

二面體的均勻原子座標 (4x4),第二個原子位於原點

dAtoms_needs_update: numpy array of bool

指出 dAtoms 是否代表 dihedraAngle

dCoordSpace: numpy array

標準化第一個多面體位置的正向和反向轉換矩陣。請參閱 dCoordSpace

dcsValid: bool

指出 dCoordSpace 是否為最新狀態

另請參閱 :meth:`build_edraArrays` 產生的屬性,以進行索引
二面體/多面體資料元素。

方法

internal_to_atom_coordinates

處理 ic 資料至殘基/原子座標;會呼叫 assemble_residues()

assemble_residues

從內部座標(平行)產生 IC_Chain 原子座標

assemble_residues_ser

從內部座標(序列)產生 IC_Residue 原子座標

atom_to_internal_coordinates

計算原子資料的二面角、角度、鍵長(內部座標)

write_SCAD

為包含鏈的內部座標資料寫入 OpenSCAD 矩陣;這是一個支援常式,請參閱 SCADIO.write_SCAD(),以產生蛋白質鏈的 OpenSCAD 描述。

distance_plot

使用選用篩選器產生原子間距離的 2D 繪圖

distance_to_internal_coordinates

從距離繪圖和二面角符號陣列計算內部座標。

make_extended

任意設定所有 psi 和 phi 骨幹角為 123 度和 -104 度。

MaxPeptideBond = 1.4

大於此值的 C-N 距離將視為鏈斷裂

ParallelAssembleResidues = True

啟用平行 internal_to_atom 演算法,對於短鏈而言速度較慢

AAsiz = 0

此鏈中的原子數(atomArray 的大小)

atomArray: array = None

AAsiz x [4] 的 float np.float64 均勻原子座標,此鏈中的所有原子。

dCoordSpace = None

[2][dihedraLen][4][4]:每個二面體的 2 個 4x4 座標空間轉換陣列。第一個 [0] 會將第一個原子轉換為標準空間,使其位於 XZ 平面上,第二個原子位於原點,第三個原子位於 +Z 軸上,第四個原子則根據二面角放置。第二個 [1] 轉換會將標準空間轉換回世界座標(PDB 檔案輸入或目前的座標)。也可以在 Dihedron 中以 .cst(正向轉換)和 .rcst(反向轉換)存取。

dcsValid = None

如果 dCoordSpace 是最新狀態則為 True。若有需要,請使用 update_dCoordSpace()

__init__(parent, verbose: bool = False) None

初始化 IC_Chain 物件,可選擇是否使用殘基/原子資料。

參數:

parent (Bio.PDB.Chain) – 此物件繼承的 Biopython Chain 物件

__deepcopy__(memo) IC_Chain

實作 IC_Chain 的深層複製。

clear_ic()

清除此鏈的殘基內部座標設定。

build_atomArray() None

從 biopython 原子建立 IC_Chain numpy 座標陣列。

另請參閱 init_edra(),以取得更完整的 IC_Chain 初始化。

輸入
self.aksetset

此鏈中的 AtomKey s

產生
self.AAsizint

鏈中的原子數 (len(akset))

self.aktupleAAsiz x AtomKeys

已排序的 akset AtomKeys

self.atomArrayIndex[AAsiz] of int

aktuple 中每個 AtomKey 的數值索引

self.atomArrayValidAAsiz x bool

如果為 True,則 atomArray 座標與內部座標一致

self.atomArrayAAsiz x np.float64[4]

均勻原子座標;執行後,Biopython Atom 座標會檢視此陣列

rak_cachedict

每個殘基的 AtomKeys 的查詢快取

build_edraArrays() None

建立鏈級的 hedra 和 dihedra 陣列。

init_edra()_hedraDict2chain() 使用。應該是私有方法,但為了文件而公開。

輸入
self.dihedraLenint

所需的 dihedra 數量

self.hedraLenint

所需的 hedra 數量

self.AAsizint

atomArray 的長度

self.hedraNdxdict

將 hedron 鍵映射到 range(hedraLen)

self.dihedraNdxdict

將 dihedron 鍵映射到 range(dihedraLen)

self.hedradict

將 Hedra 鍵映射到鏈的 Hedra

self.atomArrayAAsiz x np.float64[4]

鏈的同質原子坐標

self.atomArrayIndexdict

將 AtomKeys 映射到 atomArray

self.atomArrayValidAAsiz x bool

表示坐標是最新的

產生
self.dCoordSpace[2][dihedraLen][4][4]

轉換到/從 dihedron 坐標空間

self.dcsValiddihedraLen x bool

表示 dCoordSpace 是最新的

self.hAtomshedraLen x 3 x np.float64[4]

hCoordSpace 中的原子坐標

self.hAtomsRhedraLen x 3 x np.float64[4]

反向順序的 hAtoms(以空間換取時間)

self.hAtoms_needs_updatehedraLen x bool

表示 hAtoms, hAtoms 是最新的

self.a2h_mapAAsiz x [int …]

將 atomArrayIndex 映射到具有該原子的 hedraNdx

self.a2ha_map[hedraLen x 3]

依 hedraNdx 順序排列的 AtomNdx

self.h2aahedraLen x [int …]

將 hedraNdx 映射到 hedron 中的 atomNdx(稍後會重塑)

Hedron.ndxint

儲存在 Hedron 物件內的 self.hedraNdx 值

self.dRevdihedraLen x bool

如果為 true,則 dihedron 反向

self.dH1ndx, dH2ndx[dihedraLen]

第 1 和第 2 個 hedra 的 hedraNdx

self.h1d_maphedraLen x []

hedraNdx -> [使用 hedron 的 dihedra]

Dihedron.h1key, h2key[AtomKey …]

dihedron 的 hedron 鍵,根據需要反向

Dihedron.hedron1, hedron2Hedron

dihedron 內對 hedra 的參考

Dihedron.ndxint

Dihedron 物件內的 self.dihedraNdx 資訊

Dihedron.cst, rcstnp.float64p4][4]

Dihedron 內的 dCoordSpace 參考

self.a2da_map[dihedraLen x 4]

以 dihedraNdx 順序排列的 AtomNdx

self.d2a_map[dihedraLen x [4]]

每個二面體的 AtomNdx(重新塑形的 a2da_map)

self.dFwdbool

如果為 True,則 dihedron 未反向

self.a2d_mapAAsiz x [[dihedraNdx]

[dihedron 中原子的原子索引 0-3]],將原子索引映射到 dihedra 和其中的原子

self.dAtoms_needs_updatedihedraLen x bool

如果為 False,則 h1、h2 中的原子是最新的

assemble_residues(verbose: bool = False) None

從內部坐標產生原子坐標(向量化)。

這是 assemble_residues_ser() 的「NumPy 平行」版本。

init_atom_coords() 形成的 dihedra 開始,將每個 dihedron 從局部坐標空間轉換為蛋白質鏈坐標空間。重複直到滿足所有依賴項。

不更新 dCoordSpace,如同 assemble_residues_ser() 一樣。如果需要,請呼叫 update_dCoordSpace()。在所有原子坐標完成後,以單一操作執行會更快。

參數:

verbose (bool) – 預設為 False。報告計算已變更 dihedra 的迭代次數

產生
self.dSet: AAsiz x dihedraLen x 4

將 dihedra 中的原子映射到 atomArray

self.dSetValid[dihedraLen][4] of bool

將有效原子映射到 dihedra 以偵測 3 或 4 個有效原子

輸出坐標會寫入 atomArray。Biopython Bio.PDB.Atom 坐標是此資料的檢視。

assemble_residues_ser(verbose: bool = False, start: int | None = None, fin: int | None = None) None

從內部坐標產生 IC_Residue 原子坐標(序列)。

請參閱 assemble_residues() 以取得「numpy 平行」版本。

如果設定,則篩選 start 和 fin 之間的位置,找到每個殘基的適當開始坐標,並傳遞給 assemble()

參數:
  • verbose (bool) – 預設為 False。描述執行時間問題

  • start,fin (int) – 預設為 None。要產生坐標的子區域的開始、結束的序列位置。

init_edra(verbose: bool = False) None

建立鏈和殘基 di/hedra 結構、陣列、atomArray。

輸入

self.ordered_aa_ic_list : IC_Residue 的列表

產生
  • edra 物件,self.di/hedra(執行 _create_edra()

  • atomArray 和支援(執行 build_atomArray()

  • self.hedraLen : 結構中的 hedra 數量

  • hedraL12 : 長度、角度的 numpy 陣列(空)

  • hedraAngle ..

  • hedraL23 ..

  • self.hedraNdx : 將 hedrakeys 映射到 hedraL12 等的字典

  • self.dihedraLen : 結構中的 dihedra 數量

  • dihedraAngle ..

  • dihedraAngleRads : 角度的 np 陣列(空)

  • self.dihedraNdx : 將 dihedrakeys 映射到 dihedraAngle 的字典

init_atom_coords() None

從角度和距離設定鏈級 di/hedra 初始坐標。

初始化 hedra 和 dihedra 局部坐標空間中的原子坐標,稍後將由 dCoordSpace 矩陣適當地轉換以進行組裝。

update_dCoordSpace(workSelector: ndarray | None = None) None

計算/更新鏈 dihedra 的坐標空間轉換。

需要更新所有原子,因此會呼叫 assemble_residues()(如果所有原子都已組裝,則會立即返回)。

參數:

workSelector ([bool]) – 選取要更新的 dihedra 的選用遮罩

propagate_changes() None

追蹤 di/hedra 以使相依原子失效。

internal_to_atom_coordinates(verbose: bool = False, start: int | None = None, fin: int | None = None) None

將 IC 資料處理為殘基/原子座標。

參數:
  • verbose (bool) – 預設為 False。描述執行時間問題

  • start,fin (int) – 可選的序列位置,用於指定要處理的子區域的起始和結束位置。

注意

設定 start 或 fin 會啟用序列化的 assemble_residues_ser(),而不是(Numpy 平行化的)assemble_residues()。起始的 C-alpha 將位於原點。

atom_to_internal_coordinates(verbose: bool = False) None

計算原子資料的二面角、角度和鍵長。

產生 atomArray(透過 init_edra)、hedra 和 dihedra 的數值陣列,以及 dihedra 的座標空間轉換。

如果指定,會產生 Gly C-beta,請參閱 IC_Residue.gly_Cbeta

參數:

verbose (bool) – 預設為 False。描述執行時的問題

distance_plot(filter: ndarray | None = None) ndarray

從 atomArray 產生 2D 距離圖。

預設是計算所有原子的距離。要產生經典的 C-alpha 距離圖,請傳遞一個布林遮罩陣列,例如

atmNameNdx = internal_coords.AtomKey.fields.atm
CaSelect = [
    atomArrayIndex.get(k)
    for k in atomArrayIndex.keys()
    if k.akl[atmNameNdx] == "CA"
]
plot = cic.distance_plot(CaSelect)

或者,這將選擇所有主鏈原子

backboneSelect = [
    atomArrayIndex.get(k)
    for k in atomArrayIndex.keys()
    if k.is_backbone()
]
參數:

filter ([bool]) – 限制用於計算的原子

另請參閱

distance_to_internal_coordinates(),它需要預設的所有原子距離圖。

dihedral_signs() ndarray

取得鏈二面角陣列中每個元素的正負號陣列 (+1/-1)。

用於 distance_to_internal_coordinates()

distplot_to_dh_arrays(distplot: ndarray, dihedra_signs: ndarray) None

從距離圖載入 di/hedra 距離陣列。

從輸入的距離圖填寫 IC_Chain 陣列 hedraL12、L23、L13 和 dihedraL14 距離值陣列,並從輸入的 dihedra_signs 陣列填寫。距離圖和 di/hedra 距離陣列必須根據 IC_Chain 中的 AtomKey 對應進行索引。.hedraNdx 和 .dihedraNdx(在 IC_Chain.init_edra() 中建立)

在執行此操作之前,請呼叫 atom_to_internal_coordinates() (或至少 init_edra())以產生 a2ha_map 和 d2a_map。

已從 distance_to_internal_coordinates() 中明確移除,因此使用者可以透過其他方法填入這些鏈 di/hedra 陣列。

distance_to_internal_coordinates(resetAtoms: bool | None = True) None

從距離和手性資料計算鏈 di/hedra。

distplot_to_dh_arrays() 或其他載入器設定 hedra L12、L23、L13 和 dihedra L14 的距離屬性。

dihedraAngles 結果在最後一步乘以 dihedra_signs,以恢復在距離圖中遺失的手性資訊(結構的鏡像具有相同的距離,但二面角符號相反)。

請注意,鏈斷裂會在重建結構中造成錯誤,請使用 copy_initNCaCs() 來避免此問題

基於 Blue 的回答,在 math.stackexchange.com 上回答了「四面體的二面角與其邊長之間的關係」。另請參閱:「四面體體積的類海龍式結果」

為了完整起見,此處包含該分析中的其他值,作為註解

  • oa = hedron1 L12 如果是反向,則為 hedron1 L23

  • ob = hedron1 L23 如果是反向,則為 hedron1 L12

  • ac = hedron2 L12 如果是反向,則為 hedron2 L23

  • ab = hedron1 L13 = OA、OB 上的餘弦定律 (hedron1 L12、L23)

  • oc = hedron2 L13 = OA、AC 上的餘弦定律 (hedron2 L12、L23)

  • bc = dihedron L14

目標是 OA,沿著邊 oa 的二面角。

參數:

resetAtoms (bool) – 預設為 True。標記 di/hedra 和 atomArray 中的所有原子,以便由 internal_to_atom_coordinates() 更新。或者,將此設定為 False 並直接操作 atomArrayValiddAtoms_needs_updatehAtoms_needs_update 以減少計算量。

copy_initNCaCs(other: IC_Chain) None

從其他 IC_Chain 複製 initNCaC 原子的原子座標。

複製座標,並將初始 NCaC 和任何鏈斷裂之後的 atomArrayValid 旗標設定為 True。

如果目標有鏈斷裂,則需要用於 distance_to_internal_coordinates()(否則每個片段都會從原點開始)。

如果從其他鏈複製內部座標也很有用。

注意:IC_Residue.set_angle()IC_Residue.set_length() 會使它們的相關原子失效,因此請在呼叫此函數之前套用它們。

make_extended()

將所有 psi 和 phi 角設定為伸展構象 (123, -104)。

__annotations__ = {'atomArray': <built-in function array>, 'atomArrayIndex': "dict['AtomKey', int]", 'bpAtomArray': "list['Atom']", 'ordered_aa_ic_list': 'list[IC_Residue]'}
class Bio.PDB.internal_coords.IC_Residue(parent: Residue)

基底類別:object

此類別用於擴展 Biopython 的 Residue 物件,使其包含內部坐標資料。

參數:
parent: 此類別擴展的 biopython Residue 物件
屬性:
no_altloc: bool,預設為 False

類別變數,若為 True 則停用 ALTLOC 原子的處理,僅使用選取的原子。

accept_atoms: tuple

類別變數 accept_atoms,產生內部坐標時要使用的 PDB 原子名稱列表。預設值為:

accept_atoms = accept_mainchain + accept_hydrogens

要在內部坐標和產生的 PDB 檔案中排除氫原子,請覆寫為:

IC_Residue.accept_atoms = IC_Residue.accept_mainchain

若要僅取得主鏈原子加上醯胺質子,請使用:

IC_Residue.accept_atoms = IC_Residue.accept_mainchain + ('H',)

若要將 D 原子轉換為 H,請設定 AtomKey.d2h = True 並使用:

IC_Residue.accept_atoms = (
    accept_mainchain + accept_hydrogens + accept_deuteriums
)

請注意,accept_mainchain = accept_backbone + accept_sidechain。因此,若要產生與序列無關的構象資料,例如在二面角空間中進行結構比對,請使用:

IC_Residue.accept_atoms = accept_backbone

或設定 gly_Cbeta = True 並使用:

IC_Residue.accept_atoms = accept_backbone + ('CB',)

變更 accept_atoms 將導致 ic_rebuild 中的預設 structure_rebuild_test 失敗(顯然如此)。若有些原子被篩選掉。使用 quick=True 選項,只測試篩選過的原子坐標,以避免此問題。

目前沒有選項可以輸出使用 D 而非 H 的內部坐標。

accept_resnames: tuple

類別變數 accept_resnames,從原子產生內部坐標時,要接受的 HETATM 的三字母殘基名稱列表。HETATM 側鏈將被忽略,但正常的主鏈原子(N、CA、C、O、CB)將被包含在內。目前僅有 CYG、YCM 和 UNK;自行承擔覆寫的風險。若要產生側鏈,請將適當的項目加入 ic_data 中的 ic_data_sidechains,並在 IC_Chain.atom_to_internal_coordinates() 中提供支援。

gly_Cbeta: bool,預設為 False

類別變數 gly_Cbeta,覆寫為 True 以在 IC_Chain.atom_to_internal_coordinates() 中產生甘胺酸 CB 原子的內部坐標

IC_Residue.gly_Cbeta = True
pic_accuracy: str,預設為 “17.13f”

類別變數 pic_accuracy,設定 .pic 檔案中數值(角度、長度)的精確度。預設值設定為高,以支援重建測試中的 mmCIF 檔案精確度。如果發現重建測試因「ERROR -COORDINATES-」而失敗,且 verbose=True 顯示僅存在些微差異,請嘗試提高此值(如果僅使用 PDB 格式檔案,則可以將其降低至 9.5)。

IC_Residue.pic_accuracy = "9.5f"
residue: Biopython Residue 物件參考

此類別擴展的 Residue 物件

hedra: 由 AtomKey 的 3 元組索引的字典

構成此殘基的 Hedra

dihedra: 由 AtomKey 的 4 元組索引的字典

構成(重疊)此殘基的 Dihedra

rprev、rnext:IC_Residue 物件的列表

指向鏈中相鄰(鍵結、非遺失、可能無序)殘基的參考

atom_coords:以 AtomKey 索引的 numpy [4] 陣列字典

已移除,若需要,請使用 AtomKey 和 atomArrayIndex 來建構

ak_set:dihedra 中的 AtomKey 集合

所有與此殘基重疊的 dihedra 中的 AtomKey (請參閱 __contains__())

alt_ids: 字元列表

來自 PDB 檔案的 AltLoc ID

bfactors: 字典

從 PDB 檔案讀取的 AtomKey 索引的 B 因子

NCaCKey: AtomKey 元組列表

N、Ca、C 主鏈原子 AtomKey 元組列表;通常只有 1 個,但如果主鏈有 altloc,則會有更多個。

is20AA: bool

如果殘基是 20 種標準胺基酸之一,則為 True,根據 Residue resname 判斷

isAccept: bool

如果 is20AA 為 True 或在下方 accept_resnames 中,則為 True

rbase: tuple

殘基位置、插入碼或無、resname(如果是標準胺基酸則為 1 個字母)

cic:IC_Chain,預設為 None

父鏈 IC_Chain 物件

scale:選擇性浮點數

用於 OpenSCAD 輸出,以產生 gly_Cbeta 鍵長

方法

assemble(atomCoordsIn, resetLocation, verbose)

從內部坐標計算此殘基的原子坐標

get_angle()

傳回所傳遞鍵的角度

get_length()

傳回指定配對的鍵長

pick_angle()

尋找所傳遞鍵的 Hedron 或 Dihedron

pick_length()

尋找所傳遞 AtomKey 配對的 hedra

set_angle()

設定所傳遞鍵的角度(不更新位置)

set_length()

設定所有相關 hedra 中指定配對的鍵長

bond_rotate(delta)

依 delta 調整相關的二面角,例如旋轉 psi (N-Ca-C-N) 會以相同的量調整相鄰的 N-Ca-C-O,以避免衝突

bond_set(angle)

使用 bond_rotate 將指定的二面角設定為角度,並據此調整相關的二面角

rak(atom info)

此殘基的快取 AtomKey

accept_resnames = ('CYG', 'YCM', 'UNK')

在此處新增具有正常主鏈的非標準殘基的三字母殘基名稱。包含 CYG 是為了測試案例 4LGY(1305 個殘基的連續鏈)。可以安全地新增更多用於 N-CA-C-O 主鏈的名稱,若有更多複雜性,則需要在 accept_atomsic_data 中的 ic_data_sidechainsIC_Chain.atom_to_internal_coordinates() 中新增支援。

no_altloc: bool = False

設定為 True 可在輸入時篩選 altloc 原子,並且僅使用 Biopython 預設的 Atoms

gly_Cbeta: bool = False

在所有 Gly 殘基上建立 beta 碳原子。

將此設定為 True 會在 atom_to_internal_coordinates() 中產生 Gly C-beta 碳的內部坐標。

資料取自 2019 年 9 月 Dunbrack cullpdb_pc20_res2.2_R1.0,僅限於具有醯胺質子的結構。請參閱

PISCES:蛋白質序列篩選伺服器

「G. Wang and R. L. Dunbrack, Jr. PISCES: a protein sequence culling server. Bioinformatics, 19:1589-1591, 2003.」

從 NCACO 查詢得出的 Ala OCCACB 平均旋轉

select avg(g.rslt) as avg_rslt, stddev(g.rslt) as sd_rslt, count(*)
from
(select f.d1d, f.d2d,
(case when f.rslt > 0 then f.rslt-360.0 else f.rslt end) as rslt
from (select d1.angle as d1d, d2.angle as d2d,
(d2.angle - d1.angle) as rslt from dihedron d1,
dihedron d2 where d1.re_class='AOACACAACB' and
d2.re_class='ANACAACAO' and d1.pdb=d2.pdb and d1.chn=d2.chn
and d1.res=d2.res) as f) as g

結果

| avg_rslt          | sd_rslt          | count   |
| -122.682194862932 | 5.04403040513919 | 14098   |
pic_accuracy: str = '17.13f'
accept_backbone = ('N', 'CA', 'C', 'O', 'OXT')
accept_sidechain = ('CB', 'CG', 'CG1', 'OG1', 'OG', 'SG', 'CG2', 'CD', 'CD1', 'SD', 'OD1', 'ND1', 'CD2', 'ND2', 'CE', 'CE1', 'NE', 'OE1', 'NE1', 'CE2', 'OE2', 'NE2', 'CE3', 'CZ', 'NZ', 'CZ2', 'CZ3', 'OD2', 'OH', 'CH2', 'NH1', 'NH2')
accept_mainchain = ('N', 'CA', 'C', 'O', 'OXT', 'CB', 'CG', 'CG1', 'OG1', 'OG', 'SG', 'CG2', 'CD', 'CD1', 'SD', 'OD1', 'ND1', 'CD2', 'ND2', 'CE', 'CE1', 'NE', 'OE1', 'NE1', 'CE2', 'OE2', 'NE2', 'CE3', 'CZ', 'NZ', 'CZ2', 'CZ3', 'OD2', 'OH', 'CH2', 'NH1', 'NH2')
accept_hydrogens = ('H', 'H1', 'H2', 'H3', 'HA', 'HA2', 'HA3', 'HB', 'HB1', 'HB2', 'HB3', 'HG2', 'HG3', 'HD2', 'HD3', 'HE2', 'HE3', 'HZ1', 'HZ2', 'HZ3', 'HG11', 'HG12', 'HG13', 'HG21', 'HG22', 'HG23', 'HZ', 'HD1', 'HE1', 'HD11', 'HD12', 'HD13', 'HG', 'HG1', 'HD21', 'HD22', 'HD23', 'NH1', 'NH2', 'HE', 'HH11', 'HH12', 'HH21', 'HH22', 'HE21', 'HE22', 'HE2', 'HH', 'HH2')
accept_deuteriums = ('D', 'D1', 'D2', 'D3', 'DA', 'DA2', 'DA3', 'DB', 'DB1', 'DB2', 'DB3', 'DG2', 'DG3', 'DD2', 'DD3', 'DE2', 'DE3', 'DZ1', 'DZ2', 'DZ3', 'DG11', 'DG12', 'DG13', 'DG21', 'DG22', 'DG23', 'DZ', 'DD1', 'DE1', 'DD11', 'DD12', 'DD13', 'DG', 'DG1', 'DD21', 'DD22', 'DD23', 'ND1', 'ND2', 'DE', 'DH11', 'DH12', 'DH21', 'DH22', 'DE21', 'DE22', 'DE2', 'DH', 'DH2')
accept_atoms = ('N', 'CA', 'C', 'O', 'OXT', 'CB', 'CG', 'CG1', 'OG1', 'OG', 'SG', 'CG2', 'CD', 'CD1', 'SD', 'OD1', 'ND1', 'CD2', 'ND2', 'CE', 'CE1', 'NE', 'OE1', 'NE1', 'CE2', 'OE2', 'NE2', 'CE3', 'CZ', 'NZ', 'CZ2', 'CZ3', 'OD2', 'OH', 'CH2', 'NH1', 'NH2', 'H', 'H1', 'H2', 'H3', 'HA', 'HA2', 'HA3', 'HB', 'HB1', 'HB2', 'HB3', 'HG2', 'HG3', 'HD2', 'HD3', 'HE2', 'HE3', 'HZ1', 'HZ2', 'HZ3', 'HG11', 'HG12', 'HG13', 'HG21', 'HG22', 'HG23', 'HZ', 'HD1', 'HE1', 'HD11', 'HD12', 'HD13', 'HG', 'HG1', 'HD21', 'HD22', 'HD23', 'NH1', 'NH2', 'HE', 'HH11', 'HH12', 'HH21', 'HH22', 'HE21', 'HE22', 'HE2', 'HH', 'HH2')

變更 accept_atoms 以限制處理的原子。請參閱 IC_Residue 以了解用法。

__init__(parent: Residue) None

使用父 Biopython Residue 初始化 IC_Residue。

參數:

parent (Residue) – Biopython Residue 物件。此物件擴展的 Biopython Residue

__deepcopy__(memo)

IC_Residue 的深度複製實作。

__contains__(ak: AtomKey) bool

如果 atomkey 在此殘基中,則傳回 True。

rak(atm: str | Atom) AtomKey

快取此殘基的 AtomKey 呼叫。

__repr__() str

列印字串為父 Residue ID。

pretty_str() str

殘基 ID 的美觀字串。

set_flexible() None

對於 OpenSCAD,將 N-CA 和 CA-C 鍵標記為可撓性接頭。

請參閱 SCADIO.write_SCAD()

set_hbond() None

對於 OpenSCAD,將 H-N 和 C-O 鍵標記為氫鍵(磁鐵)。

請參閱 SCADIO.write_SCAD()

clear_transforms()

在 assemble() 之前,使二面角座標空間屬性失效。

座標空間屬性為 Dihedron.cst 和 .rcst,以及 IC_Chain.dCoordSpace

assemble(resetLocation: bool = False, verbose: bool = False) dict[AtomKey, array] | dict[tuple[AtomKey, AtomKey, AtomKey], array] | None

從內部座標計算此殘基的原子座標。

這是 assemble_residues_ser() 序列版本的 IC_Residue 部分,請參閱 assemble_residues(),以了解在 IC_Chain 層級上運作的 numpy 向量化方法。

從 N-CA-C 和 N-CA-CB 四面體開始,連接準備好的二面角,計算主鏈和側鏈原子的蛋白質空間座標

在每個二面角上設定正向和反向轉換,將前三個原子的蛋白質座標轉換為二面角空間座標(請參閱 IC_Chain.dCoordSpace

呼叫 init_atom_coords() 以更新在此之前修改的任何二/四面體,這裡僅將二面體組裝到蛋白質座標空間中。

演算法

形成雙端佇列,從 c-ca-n、o-c-ca、n-ca-cb、n-ca-c 開始。

如果 resetLocation=True,則使用生成二面體的初始座標作為 n-ca-c 的初始位置(在二面角座標空間中)

當佇列不為空時

取得 3 個原子的四面體鍵

對於每個以四面體鍵(二面角的第 1 個四面體)開頭的二面角

如果所有 4 個原子都已有座標

將第 2 個四面體鍵添加到佇列末尾

否則,如果前 3 個原子有座標

計算正向和反向轉換,將前 3 個原子從/到二面角的初始座標空間

使用反向轉換從二面角的初始座標取得第 4 個原子在目前座標中的位置

將第 2 個四面體鍵添加到佇列末尾

否則

排序失敗,將四面體鍵放回佇列末尾,並希望下次我們有前 3 個原子位置(不應該發生)

迴圈終止(佇列排空),因為沒有開始任何二面角的四面體鍵在沒有動作的情況下被移除

參數:

resetLocation (bool) – 預設為 False。- 忽略開始位置並定向的選項,使初始 N-Ca-C 四面體位於原點。

傳回:

AtomKey -> 相對於前一個殘基的蛋白質空間中殘基的均勻原子座標的字典

此外,如同 assemble_residues() 一樣,直接更新 IC_Chain.atomArray

split_akl(lst: tuple[AtomKey, ...] | list[AtomKey], missingOK: bool = False) list[tuple[AtomKey, ...]]

取得此殘基的 AtomKeys (ak_set),用於 AtomKeys 的通用清單。

更改和/或擴展「通用」AtomKeys 的清單(例如「N、C、C」),使其特定於此殘基的 altloc 等,例如「(N-Ca_A_0.3-C, N-Ca_B_0.7-C)」

給定四面體或二面角的 AtomKeys 清單,
傳回

在此殘基中具有 id3_dh 的匹配原子鍵清單(如果佔用率 != 1.00,則 ak 可能會變更)

為所有原子 altloc 擴展的匹配原子鍵的多個清單

如果任何 atom_coord(ak) 遺失且未 missingOK,則為空清單

參數:
  • lst (list) – AtomKeys 的清單 [3] 或 [4]。要與此殘基的特定 AtomKeys 匹配的非 altloc AtomKeys

  • missingOK (bool) – 預設為 False,請參閱上方。

atom_sernum = None
atom_chain = None
pdb_residue_string() str

以字串形式產生此殘基的 PDB ATOM 記錄。

PDBIO.py 中未公開的功能的便利方法。如果不是 None,則遞增 IC_Residue.atom_sernum

參數:
  • IC_Residue.atom_sernum – 類別變數,預設為 None。如果不是 None,則覆寫並遞增原子序號

  • IC_Residue.atom_chain – 類別變數。如果不是 None,則覆寫原子鏈 ID

待辦事項

移動到 PDBIO

pic_flags = (1, 2, 4, 8, 16, 32, 64, 128, 256, 512, 496, 525, 1021, 1024, 2048, 4096, 7168, 8192, 16384)

PICIO.write_PIC() 使用,以控制要預設的數值類別。

picFlagsDefault = 31744

預設為所有二面角 + 初始 tau 原子 + bFactors。

picFlagsDict = {'all': 7168, 'bFactors': 16384, 'chi': 496, 'chi1': 16, 'chi2': 32, 'chi3': 64, 'chi4': 128, 'chi5': 256, 'classic': 1021, 'classic_b': 525, 'hedra': 1024, 'initAtoms': 8192, 'omg': 2, 'phi': 4, 'pomg': 512, 'primary': 2048, 'psi': 1, 'secondary': 4096, 'tau': 8}

根據需要使用的 pic_flags 值的字典。

pick_angle(angle_key: tuple[AtomKey, AtomKey, AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey] | str) Hedron | Dihedron | None

取得 angle_key 的四面體或二面角。

參數:

angle_key

  • 3 或 4 個 AtomKeys 的元組

  • 由 : 分隔的原子名稱字串 ('CA')

  • 由 ':' 分隔的 [-1, 0, 1]<原子名稱> 字串。-1 是上一個殘基,0 是此殘基,1 是下一個殘基

  • psi、phi、omg、omega、chi1、chi2、chi3、chi4、chi5

  • tau (N-CA-C 角),請參閱 Richardson1981

  • AtomKeys 的元組僅適用於存取替代無序原子

請注意,殘基的 phi 和 omega 二面角,以及構成它們的四面體(包括 N:Ca:C tau 四面體)都儲存在 n-1 個二/四面體集合中;此重疊在此處處理,但如果直接存取可能會出現問題。

下列 print 指令是等效的(對於 chi2 具有非碳原子的側鏈除外)

ric = r.internal_coord
print(
    r,
    ric.get_angle("psi"),
    ric.get_angle("phi"),
    ric.get_angle("omg"),
    ric.get_angle("tau"),
    ric.get_angle("chi2"),
)
print(
    r,
    ric.get_angle("N:CA:C:1N"),
    ric.get_angle("-1C:N:CA:C"),
    ric.get_angle("-1CA:-1C:N:CA"),
    ric.get_angle("N:CA:C"),
    ric.get_angle("CA:CB:CG:CD"),
)

請參閱 ic_data.py,以了解列舉側鏈角度和不跨越肽鍵的主鏈角度中的原子的詳細資訊。使用 's' 表示目前殘基 ('self'),'n' 表示下一個殘基,跨越(重疊)角度為

(sN, sCA, sC, nN)   # psi
(sCA, sC, nN, nCA)  # omega i+1
(sC, nN, nCA, nC)   # phi i+1
(sCA, sC, nN)
(sC, nN, nCA)
(nN, nCA, nC)       # tau i+1
傳回:

匹配的四面體、二面角或 None。

get_angle(angle_key: tuple[AtomKey, AtomKey, AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey] | str) float | None

取得指定鍵的二面角或多面角角度。

關於鍵的詳細規範,請參閱 pick_angle()

set_angle(angle_key: tuple[AtomKey, AtomKey, AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey] | str, v: float, overlap=True)

設定指定鍵的二面角或多面角角度。

如果角度是 Dihedronoverlap 為 True (預設值),則重疊的二面角也會做相應的變更。重疊是 ic_data_create_edra() 中蛋白質鏈定義的結果(例如,psi 與 N-CA-C-O 重疊)。

對於以下情況,預設值 overlap=True 可能就是您想要的:set_angle(“chi1”, val)

當處理鏈或殘基中的所有二面角(例如從另一個結構複製)時,預設值可能不是您想要的,因為重疊的二面角也可能在該集合中。

注意:設定例如 PRO chi2 是允許的,不會產生錯誤或警告!

關於 angle_key 的詳細規範,請參閱 pick_angle()。若要將二面角變更一定的角度,請參閱 bond_rotate()

參數:
  • angle_key – 角度識別碼。

  • v (float) – 以度為單位的新角度(結果調整為 +/-180)。

  • overlap (bool) – 預設值為 True。根據需要修改重疊的二面角。

bond_rotate(angle_key: tuple[AtomKey, AtomKey, AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey] | str, delta: float)

將一組重疊的二面角旋轉 delta 度。

將二面角角度變更給定的 delta 值,即 new_angle = current_angle + delta。數值會經過調整,使 new_angle 在 +/-180 範圍內。

變更重疊的二面角,如同 set_angle()

關於鍵的詳細規範,請參閱 pick_angle()

bond_set(angle_key: tuple[AtomKey, AtomKey, AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey] | str, val: float)

將二面角設定為 val,並將重疊的二面角更新相同的量。

set_angle() 類似,保留此方法是為了相容性。與 set_angle() 不同的是,此方法僅適用於二面角,且沒有不更新重疊二面角的選項。

關於鍵的詳細規範,請參閱 pick_angle()

pick_length(ak_spec: str | tuple[AtomKey, AtomKey]) tuple[list[Hedron] | None, tuple[AtomKey, AtomKey] | None]

取得包含指定原子對的多面體列表。

參數:

ak_spec

  • 兩個 AtomKey 的元組

  • 字串:兩個原子名稱以 ‘:’ 分隔,例如 ‘N:CA’,並可選擇指定相對於自身的定位符,例如 ‘-1C:N’ 代表前一個胜肽鍵。定位符為 -1、0、1。

以下是等效的

ric = r.internal_coord
print(
    r,
    ric.get_length("0C:1N"),
)
print(
    r,
    None
    if not ric.rnext
    else ric.get_length((ric.rak("C"), ric.rnext[0].rak("N"))),
)

如果在目前殘基上找不到原子,則會在 rprev[0] 中尋找,以處理類似 Gly N:CA 的情況。如需更精細的控制,請直接存取 IC_Chain.hedra

傳回:

包含指定原子對的晶體結構列表,以 AtomKeys 元組形式表示

get_length(ak_spec: str | tuple[AtomKey, AtomKey]) float | None

取得指定原子對的鍵長。

請參閱 pick_length() 以了解 ak_spec 和詳細資訊。

set_length(ak_spec: str | tuple[AtomKey, AtomKey], val: float) None

設定指定原子對的鍵長。

請參閱 pick_length() 以了解 ak_spec。

applyMtx(mtx: array) None

將矩陣套用至此 IC_Residue 的 atom_coords。

__annotations__ = {'_AllBonds': <class 'bool'>, 'ak_set': 'set[AtomKey]', 'akc': 'dict[Union[str, Atom], AtomKey]', 'alt_ids': 'Union[list[str], None]', 'bfactors': 'dict[str, float]', 'cic': 'IC_Chain', 'dihedra': 'dict[DKT, Dihedron]', 'gly_Cbeta': <class 'bool'>, 'hedra': 'dict[HKT, Hedron]', 'no_altloc': <class 'bool'>, 'pic_accuracy': <class 'str'>, 'rnext': 'list[IC_Residue]', 'rprev': 'list[IC_Residue]'}
class Bio.PDB.internal_coords.Edron(*args: list[AtomKey] | tuple[AtomKey, AtomKey, AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey], **kwargs: str)

基底類別:object

Hedron 和 Dihedron 類別的基底類別。

支援基於 AtomKeys 列表的豐富比較。

屬性:
atomkeys:元組

定義此二面體/多面體的 3 個(多面體)或 4 個(二面體)AtomKey

id:字串

此二面體/多面體的 AtomKeys 的 ‘:’ 連接字串

needs_update:布林值

表示二面體/多面體局部 atom_coords 並未反映多面體局部座標空間中目前的二面體/多面體角度和長度值

e_class:字串

組成二面體/多面體的原子序列(無位置或殘基),用於統計

re_class:字串

組成二面體/多面體的殘基、原子序列,用於統計

cre_class:字串

組成二面體/多面體的共價半徑類別序列,用於統計

edron_re:編譯的 regex(類別屬性)

一個編譯的正規表示式,用於比對 Hedron 和 Dihedron 物件的字串 ID

cic:IC_Chain 參考

包含此多面體的鏈內部座標物件

ndx:整數

進入 IC_Chain 層級 numpy 資料陣列的二面體/多面體索引。在 IC_Chain.init_edra() 中設定

rc:整數

此 edron 中涉及的殘基數量

方法

gen_key([AtomKey, …] 或 AtomKey, …) (靜態方法)

產生 AtomKey ID 的 ‘:’ 連接字串

is_backbone()

如果所有 atomkeys 原子都是 N、Ca、C 或 O,則傳回 True

edron_re = re.compile('^(?P<pdbid>\\w+)?\\s(?P<chn>[\\w|\\s])?\\s(?P<a1>[\\w\\-\\.]+):(?P<a2>[\\w\\-\\.]+):(?P<a3>[\\w\\-\\.]+)(:(?P<a4>[\\w\\-\\.]+))?\\s+(((?P<len12>\\S+)\\s+(?P<angle>\\S+)\\s+(?P<len23>\\S+)\\s*$)|((?P<)

一個編譯的正規表示式,用於比對 Hedron 和 Dihedron 物件的字串 ID

static gen_key(lst: list[AtomKey]) str

從輸入產生 ‘:’ 連接的 AtomKey 字串。

從 (2_A_C, 3_P_N, 3_P_CA) 產生 ‘2_A_C:3_P_N:3_P_CA’ :param list lst: AtomKey 物件的列表

static gen_tuple(akstr: str) tuple

為 ‘:’ 連接的 AtomKey 字串產生 AtomKey 元組。

從 ‘2_A_C:3_P_N:3_P_CA’ 產生 (2_A_C, 3_P_N, 3_P_CA) :param str akstr: ‘:’ 分隔的 AtomKey 字串

__init__(*args: list[AtomKey] | tuple[AtomKey, AtomKey, AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey], **kwargs: str) None

使用 AtomKeys 序列初始化 Edron。

可接受的輸入

[ AtomKey, … ]:AtomKeys 的列表 AtomKey, …:AtomKeys 序列作為 args {‘a1’: str, ‘a2’: str, … }:AtomKeys 的字典,如 ‘a1’、‘a2’ …

__deepcopy__(memo)

Edron 的深層複製實作。

__contains__(ak: AtomKey) bool

如果 atomkey 在此 edron 中,則傳回 True。

is_backbone() bool

僅包含 N、C、CA、O、H 原子時回報 True。

__repr__() str

AtomKeys 的元組是預設的 repr 字串。

__hash__() int

雜湊值在初始化時從 atomkeys 元組計算得出。

__eq__(other: object) bool

測試是否相等。

__ne__(other: object) bool

測試是否不相等。

__gt__(other: object) bool

測試是否大於。

__ge__(other: object) bool

測試是否大於或等於。

__lt__(other: object) bool

測試是否小於。

__le__(other: object) bool

測試是否小於或等於。

class Bio.PDB.internal_coords.Hedron(*args: list[AtomKey] | tuple[AtomKey, AtomKey, AtomKey], **kwargs: str)

基礎類別: Edron

表示形成平面的三個相連原子的類別。

包含局部座標空間中的原子座標:中心原子位於原點,一個終端原子位於 XZ 平面上,另一個終端原子位於 +Z 軸上。儲存為兩個方向,其中第 3 個(正向)或第 1 個(反向)原子位於 +Z 軸上。有關正向和反向方向的使用,請參閱 Dihedron

屬性:
len12: float

第一個和第二個原子之間的距離

len23: float

第二個和第三個原子之間的距離

angle: float

hedron 中三個原子形成的角度(度)

xrh_class: string

僅適用於跨越 2 個殘基的 hedron,對於僅貢獻一個原子的殘基,將具有 'X'

方法

get_length()

取得指定原子對的鍵長

set_length()

設定指定原子對的鍵長

angle(), len12(), len23()

相關屬性的 setter(角度以度為單位)

__init__(*args: list[AtomKey] | tuple[AtomKey, AtomKey, AtomKey], **kwargs: str) None

使用 AtomKeys 序列和 kwargs 初始化 Hedron。

可接受的輸入

如同 Edron,加上可選的 'len12'、'angle'、'len23' 關鍵字值。

__repr__() str

列印 Hedron 物件的字串。

property angle: float

取得此 hedron 的角度。

property len12

取得 Hedron 的第一個長度。

property len23: float

取得 Hedron 的第二個長度。

get_length(ak_tpl: tuple[AtomKey, AtomKey]) float | None

取得指定原子對的鍵長。

參數:

ak_tpl (tuple) – AtomKeys 的元組。此 Hedron 中的原子對

set_length(ak_tpl: tuple[AtomKey, AtomKey], newLength: float)

設定指定原子對的鍵長;設定 needs_update。

參數:

.ak_tpl (tuple) – AtomKeys 的元組。此 Hedron 中的原子對

class Bio.PDB.internal_coords.Dihedron(*args: list[AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey], **kwargs: str)

基礎類別: Edron

表示形成二面角的四個相連原子的類別。

屬性:
angle: float

以度為單位測量或指定二面角;首選使用 IC_Residue.bond_set() 進行設定

hedron1, hedron2:Hedron 物件參考

形成二面角的兩個 hedra

h1key, h2key:AtomKeys 的元組

hedron1 和 hedron2 的雜湊鍵

id3,id32:AtomKeys 的元組

構成二面角的前 3 個和後 3 個原子;hxkey 的順序可能不同

ric:IC_Residue 物件參考

包含此二面角的 IC_Residue 物件

reverse:bool

表示二面角中原子的順序與 hedra 中原子的順序相反

primary:bool

如果這是 psi、phi、omega 或側鏈 chi 角,則為 True

pclass:字串(主要角度類別)

根據命名法(psi、omega、phi),以 X 表示相鄰殘基的 re_class

cst, rcst:numpy [4][4] 陣列

轉換到(cst)和從(rcst)二面角坐標空間,該空間以原子 2(Hedron 1 中心原子)為原點定義。請參閱 IC_Chain.dCoordSpace

方法

angle()

用於取得/設定二面角角度(以度為單位);建議使用 IC_Residue.bond_set()

bits()

返回 IC_Residue.pic_flags 位元遮罩,用於二面角 psi、omega 等

__init__(*args: list[AtomKey] | tuple[AtomKey, AtomKey, AtomKey, AtomKey], **kwargs: str) None

使用 AtomKeys 序列和可選的二面角初始化 Dihedron。

可接受的輸入

與 Edron 相同,加上可選的「dihedral」關鍵字角度值。

__repr__() str

列印 Dihedron 物件的字串。

property angle: float

取得二面角。

static angle_dif(a1: float | ndarray, a2: float | ndarray)

取得兩個 +/- 180 度角之間的角度差。

https://stackoverflow.com/a/36001014/2783487

static angle_avg(alst: list, in_rads: bool = False, out_rads: bool = False)

取得 +/-180 度角列表的平均值。

參數:
  • alst (List) – 要平均的角度列表

  • in_rads (bool) – 輸入值以弧度表示

  • out_rads (bool) – 以弧度回報結果

static angle_pop_sd(alst: list, avg: float)

取得 +/-180 度角列表的母體標準差。

應該是樣本標準差,但避免 len(alst)=1 -> 除以 0

difference(other: Dihedron) float

取得此角度和另一個 +/- 180 度角之間的角度差。

bits() int

取得 IC_Residue.pic_flags 位元遮罩,用於自身是否為 psi、omg、phi、pomg、chiX。

class Bio.PDB.internal_coords.AtomKey(*args: IC_Residue | Atom | list | dict | str, **kwargs: str)

基底類別:object

用於參考原子坐標的字典鍵的類別。

AtomKeys 會一起擷取殘基和無序資訊,並為 .pic 檔案提供不含空格的字串鍵。

支援豐富的比較和多種實例化的方式。

AtomKeys 包含

殘基位置 (respos)、插入碼 (icode)、1 或 3 個字元的殘基名稱 (resname)、原子名稱 (atm)、altloc (altloc) 和佔有率 (occ)

使用 AtomKey.fields 按名稱取得感興趣的組件的索引

使用 AtomKeys 從 IC_Chain atomArray 和 atomArrayIndex 取得 C-alpha 原子

atmNameNdx = internal_coords.AtomKey.fields.atm
CaSelection = [
    atomArrayIndex.get(k)
    for k in atomArrayIndex.keys()
    if k.akl[atmNameNdx] == "CA"
]
AtomArrayCa = atomArray[CaSelection]

取得鏈中的所有苯丙氨酸原子

resNameNdx = internal_coords.AtomKey.fields.resname
PheSelection = [
    atomArrayIndex.get(k)
    for k in atomArrayIndex.keys()
    if k.akl[resNameNdx] == "F"
]
AtomArrayPhe = atomArray[PheSelection]

如果為 20 個標準殘基之一,「resname」將為大寫的 1 個字母氨基酸代碼,否則為提供的 3 個字母代碼。作為輸入提供,或從 IC_Residue 的 .rbase 屬性讀取。

屬性:
akl: tuple

AtomKey 的所有六個欄位

fieldNames: tuple (類別屬性)

索引鍵位置到名稱的對應

fields: namedtuple (類別屬性)

欄位名稱到索引位置的對應。

id:字串

以「_」連接的 AtomKey 欄位,不包含「None」欄位

atom_re: 編譯的正規表示式 (類別屬性)

一個編譯的正規表示式,可比對金鑰的字串形式

d2h: bool (類別屬性) 預設為 False

如果為 True,則在輸入時將 D 原子轉換為 H;還必須修改 IC_Residue.accept_atoms

missing: bool 預設為 False

從字串 __init__ 的 AtomKey 可能遺失,請設定此旗標以註記此問題。IC_Residue.rak() 設定

ric: IC_Residue 預設為 None

如果使用 IC_Residue 初始化,則此值會參考 IC_residue

方法

altloc_match(other)

如果此 AtomKey 符合其他 AtomKey,但不包含 altloc 和佔有率欄位,則傳回 True

is_backbone()

如果原子為 N、CA、C、O 或 H,則傳回 True

atm()

傳回原子名稱,例如 N、CA、CB 等。

cr_class()

傳回共價半徑類別,例如 Csb

atom_re = re.compile('^(?P<respos>-?\\d+)(?P<icode>[A-Za-z])?_(?P<resname>[a-zA-Z]+)_(?P<atm>[A-Za-z0-9]+)(?:_(?P<altloc>\\w))?(?:_(?P<occ>-?\\d\\.\\d+?))?$')

預先編譯的正規表示式,可比對 AtomKey 字串。

fieldNames = ('respos', 'icode', 'resname', 'atm', 'altloc', 'occ')
fields = (0, 1, 2, 3, 4, 5)

使用此 namedtuple 來存取 AtomKey 欄位。請參閱 AtomKey

d2h = False

設定為 True 以在輸入時將 D 氘轉換為 H 氫。

__init__(*args: IC_Residue | Atom | list | dict | str, **kwargs: str) None

使用殘基和原子數據初始化 AtomKey。

可接受輸入的範例

(<IC_Residue>, 'CA', ...)    : IC_Residue with atom info
(<IC_Residue>, <Atom>)       : IC_Residue with Biopython Atom
([52, None, 'G', 'CA', ...])  : list of ordered data fields
(52, None, 'G', 'CA', ...)    : multiple ordered arguments
({respos: 52, icode: None, atm: 'CA', ...}) : dict with fieldNames
(respos: 52, icode: None, atm: 'CA', ...) : kwargs with fieldNames
52_G_CA, 52B_G_CA, 52_G_CA_0.33, 52_G_CA_B_0.33  : id strings
__deepcopy__(memo)

AtomKey 的深層複製實作。

__repr__() str

來自 id 的 Repr 字串。

__hash__() int

從 akl 元組在初始化時計算的雜湊值。

altloc_match(other: AtomKey) bool

測試 AtomKey 與其他鍵是否匹配,忽略佔有率和替選位置。

is_backbone() bool

如果是 N、C、CA、O 或 H,則返回 True。

atm() str

返回原子名稱:N、CA、CB、O 等。

cr_class() str | None

返回原子的共價半徑類別,如果沒有則返回 None。

__ne__(other: object) bool

測試是否不相等。

__eq__(other: object) bool

測試是否相等。

__gt__(other: object) bool

測試是否大於。

__ge__(other: object) bool

測試是否大於或等於。

__lt__(other: object) bool

測試是否小於。

__le__(other: object) bool

測試是否小於或等於。

Bio.PDB.internal_coords.set_accuracy_95(num: float) float

將浮點數精確度降低至 9.5 (xxxx.xxxxx)。

IC_Residue 類別在寫入 PIC 和 SCAD 檔案時使用。

參數:

num (float) – 輸入數字

傳回:

具有指定精確度的浮點數

exception Bio.PDB.internal_coords.HedronMatchError

基礎:Exception

無法在殘基中找到給定鍵的四面體。

exception Bio.PDB.internal_coords.MissingAtomError

基礎:Exception

四面體或二面體缺少原子坐標。