Gemmi C++ API
Loading...
Searching...
No Matches
modify.hpp
Go to the documentation of this file.
1// Copyright 2017-2021 Global Phasing Ltd.
2//
3// Modify various properties of the model.
4
5// For modifications that depend on entities or connectivity see polyheur.hpp.
6
7#ifndef GEMMI_MODIFY_HPP_
8#define GEMMI_MODIFY_HPP_
9
10#include "model.hpp"
11#include "util.hpp" // for vector_remove_if
12#include <set>
13
14namespace gemmi {
15
17template<class T> void remove_alternative_conformations(T& obj) {
18 for (auto& child : obj.children())
20}
21template<> inline void remove_alternative_conformations(Chain& chain) {
22 std::set<SeqId> seqids;
23 for (size_t i = 0; i < chain.residues.size(); ) {
24 if (seqids.insert(chain.residues[i].seqid).second)
25 ++i;
26 else
27 chain.residues.erase(chain.residues.begin() + i);
28 }
29 for (Residue& residue : chain.residues) {
30 std::set<std::string> names;
31 for (size_t i = 0; i < residue.atoms.size(); ) {
32 Atom& atom = residue.atoms[i];
33 atom.altloc = '\0';
34 if (names.insert(atom.name).second)
35 ++i;
36 else
37 residue.atoms.erase(residue.atoms.begin() + i);
38 }
39 }
40}
41
43template<class T> void remove_hydrogens(T& obj) {
44 for (auto& child : obj.children())
46}
47template<> inline void remove_hydrogens(Residue& res) {
48 vector_remove_if(res.atoms, [](const Atom& a) {
49 return a.element == El::H || a.element == El::D;
50 });
51}
52
56template<class T> void assign_b_iso(T& obj, float b_min, float b_max) {
57 for (auto& child : obj.children())
59}
60template<> inline void assign_b_iso(Atom& atom, float b_min, float b_max) {
61 atom.b_iso = std::min(std::max(atom.b_iso, b_min), b_max);
62}
63
65template<class T> void remove_anisou(T& obj) {
66 for (auto& child : obj.children())
68}
69template<> inline void remove_anisou(Atom& atom) {
70 atom.aniso = {0, 0, 0, 0, 0, 0};
71}
72
74template<class T> void ensure_anisou(T& obj) {
75 for (auto& child : obj.children())
77}
78template<> inline void ensure_anisou(Atom& atom) {
79 if (!atom.aniso.nonzero()) {
80 float u = float(1. / gemmi::u_to_b() * atom.b_iso);
81 atom.aniso = {u, u, u, 0.f, 0.f, 0.f};
82 }
83}
84
86template<class T> void transform_pos_and_adp(T& obj, const Transform& tr) {
87 for (auto& child : obj.children())
89}
90template<> inline void transform_pos_and_adp(Atom& atom, const Transform& tr) {
91 atom.pos = Position(tr.apply(atom.pos));
92 if (atom.aniso.nonzero())
93 atom.aniso = atom.aniso.transformed_by<float>(tr.mat);
94}
95
97inline void assign_serial_numbers(Model& model, bool numbered_ter=false) {
98 int serial = 0;
99 for (Chain& chain : model.chains)
100 for (Residue& res : chain.residues) {
101 for (Atom& atom : res.atoms)
102 atom.serial = ++serial;
103 if (numbered_ter && res.entity_type == EntityType::Polymer &&
104 (&res == &chain.residues.back() || (&res + 1)->entity_type != EntityType::Polymer))
105 ++serial;
106 }
107}
108inline void assign_serial_numbers(Structure& st, bool numbered_ter=false) {
109 for (Model& model : st.models)
110 assign_serial_numbers(model, numbered_ter);
111}
112
113
114template<typename Func>
116 for (Connection& con : st.connections) {
117 func(con.partner1);
118 func(con.partner2);
119 }
120 for (CisPep& cispep : st.cispeps) {
121 func(cispep.partner_c);
122 func(cispep.partner_n);
123 }
124 for (Helix& helix : st.helices) {
125 func(helix.start);
126 func(helix.end);
127 }
128 for (Sheet& sheet : st.sheets)
129 for (Sheet::Strand& strand : sheet.strands) {
130 func(strand.start);
131 func(strand.end);
132 func(strand.hbond_atom2);
133 func(strand.hbond_atom1);
134 }
135}
136
137
138inline void rename_chain(Structure& st, const std::string& old_name,
139 const std::string& new_name) {
141 if (aa.chain_name == old_name)
142 aa.chain_name = new_name;
143 });
144 for (RefinementInfo& ri : st.meta.refinement)
145 for (TlsGroup& tls : ri.tls_groups)
146 for (TlsGroup::Selection& sel : tls.selections)
147 if (sel.chain == old_name)
148 sel.chain = new_name;
149 for (Model& model : st.models)
150 for (Chain& chain : model.chains)
151 if (chain.name == old_name)
152 chain.name = new_name;
153}
154
155inline void rename_atom_names(Structure& st, const std::string& res_name,
156 const std::map<std::string, std::string>& old_new) {
157 auto update = [&old_new](std::string& name) {
158 auto it = old_new.find(name);
159 if (it != old_new.end())
160 name = it->second;
161 };
163 if (aa.res_id.name == res_name)
164 update(aa.atom_name);
165 });
166 for (Model& model : st.models)
167 for (Chain& chain : model.chains)
168 for (Residue& res : chain.residues)
169 if (res.name == res_name) {
170 for (Atom& atom : res.atoms)
171 update(atom.name);
172 }
173}
174
175
177 for (size_t i = res.atoms.size(); i-- != 0; ) {
178 Atom& atom = res.atoms[i];
179 float d_fraction = atom.fraction;
180 if (atom.element == El::H && d_fraction > 0) {
181 if (d_fraction >= 1) {
182 atom.element = El::D;
183 if (atom.name[0] == 'H')
184 atom.name[0] = 'D';
185 } else {
186 int alt_offset = atom.altloc;
187 if (alt_offset) {
188 alt_offset -= 'A';
189 // we don't expect 4+ altlocs - ignore such cases
191 continue;
192 }
193 atom.altloc = 'A' + alt_offset;
194 float d_occ = atom.occ * d_fraction;
195 atom.occ *= (1 - d_fraction);
196 auto deut = res.atoms.insert(res.atoms.begin() + i + 1, atom);
197 deut->altloc = 'D' + alt_offset;
198 deut->element = El::D;
199 deut->occ = d_occ;
200 if (deut->name[0] == 'H')
201 deut->name[0] = 'D';
202 }
203 }
204 }
205}
206
208 bool found = false;
209 for (auto d = res.atoms.end(); d-- != res.atoms.begin(); )
210 if (d->element == El::D) {
211 found = true;
212 auto h = res.atoms.begin();
213 for (; h != res.atoms.end(); ++h)
214 if (h->element == El::H && h->pos.approx(d->pos, 1e-9))
215 break;
216 if (h != res.atoms.end()) {
217 h->occ += d->occ;
218 h->fraction = h->occ > 0.f ? d->occ / h->occ : 0.f;
219 if (h->altloc) {
220 bool keep_altloc = false;
221 for (auto i = res.atoms.begin(); i != res.atoms.end(); ++i)
222 if (i != d && i != h && (i->name == h->name || i->name == d->name))
223 keep_altloc = true;
224 if (!keep_altloc)
225 h->altloc = '\0';
226 }
227 res.atoms.erase(d);
228 } else {
229 d->element = El::H;
230 d->fraction = 1;
231 // Atom name is left unchanged. prepare_topology() first calls this
232 // function and then conditionally changes the name (Dxx -> Hxx).
233 }
234 }
235 return found;
236}
237
244 if (st.has_d_fraction == store_fraction)
245 return;
246 st.has_d_fraction = false;
247 for (Model& model : st.models)
248 for (Chain& chain : model.chains)
249 for (Residue& res : chain.residues)
250 if (store_fraction) {
252 st.has_d_fraction = true;
253 } else {
255 }
256}
257
260 if (!st.cell.explicit_matrices || !st.cell.is_crystal())
261 return;
262 Transform orig_frac = st.cell.frac;
263 st.cell.explicit_matrices = false;
264 st.cell.calculate_properties();
265 Transform tr = st.cell.orth.combine(orig_frac);
266 Transform tr_inv = tr.inverse();
267 st.has_origx = true;
268 st.origx = tr_inv.combine(st.origx);
269 for (NcsOp& ncsop : st.ncs)
270 ncsop.tr = tr.combine(ncsop.tr).combine(tr_inv);
272}
273
274} // namespace gemmi
275#endif
void assign_serial_numbers(Model &model, bool numbered_ter=false)
set atom site serial numbers to 1, 2, ..., optionally leaving gaps for TERs
Definition modify.hpp:97
void remove_anisou(T &obj)
Remove anisotropic ADP.
Definition modify.hpp:65
void remove_alternative_conformations(T &obj)
Remove alternative conformations.
Definition modify.hpp:17
void replace_d_fraction_with_altlocs(Residue &res)
Definition modify.hpp:176
void vector_remove_if(std::vector< T > &v, F &&condition)
Definition util.hpp:266
bool replace_deuterium_with_fraction(Residue &res)
Definition modify.hpp:207
void process_addresses(Structure &st, Func func)
Definition modify.hpp:115
void rename_chain(Structure &st, const std::string &old_name, const std::string &new_name)
Definition modify.hpp:138
void transform_pos_and_adp(T &obj, const Transform &tr)
apply Transform to both atom's position and ADP
Definition modify.hpp:86
void remove_hydrogens(T &obj)
Remove hydrogens.
Definition modify.hpp:43
void store_deuterium_as_fraction(Structure &st, bool store_fraction)
Hydrogens modelled as H/D mixture (altlocs H and D with the same position and ADP,...
Definition modify.hpp:243
void rename_atom_names(Structure &st, const std::string &res_name, const std::map< std::string, std::string > &old_new)
Definition modify.hpp:155
void standardize_crystal_frame(Structure &st)
Convert coordinates to the standard coordinate system for the unit cell.
Definition modify.hpp:259
void assign_b_iso(T &obj, float b_min, float b_max)
Set isotropic ADP to the range (b_min, b_max).
Definition modify.hpp:56
void ensure_anisou(T &obj)
Set absent ANISOU to value from B_iso.
Definition modify.hpp:74
constexpr double u_to_b()
Definition math.hpp:29
Represents atom site in macromolecular structure (~100 bytes).
Definition model.hpp:112
SMat33< float > aniso
Definition model.hpp:127
char altloc
Definition model.hpp:115
float b_iso
Definition model.hpp:126
Element element
Definition model.hpp:117
float fraction
Definition model.hpp:122
float occ
Definition model.hpp:124
Position pos
Definition model.hpp:123
std::string name
Definition model.hpp:114
std::vector< Residue > residues
Definition model.hpp:465
std::vector< Chain > chains
Definition model.hpp:700
Non-crystallographic symmetry operation (such as in the MTRIXn record)
Definition unitcell.hpp:120
Coordinates in Angstroms - orthogonal (Cartesian) coordinates.
Definition unitcell.hpp:32
std::vector< Atom > atoms
Definition model.hpp:188
bool nonzero() const
Definition math.hpp:273
SMat33< Real > transformed_by(const Mat33 &m) const
Definition math.hpp:319
Vec3 apply(const Vec3 &x) const
Definition math.hpp:380
Transform combine(const Transform &b) const
Definition math.hpp:382
Transform inverse() const
Definition math.hpp:375