7#ifndef GEMMI_MODIFY_HPP_
8#define GEMMI_MODIFY_HPP_
23 for (
size_t i = 0;
i < chain.
residues.size(); ) {
30 std::set<std::string>
names;
31 for (
size_t i = 0;
i < residue.atoms.size(); ) {
32 Atom& atom = residue.atoms[
i];
37 residue.atoms.erase(residue.atoms.begin() +
i);
49 return a.element == El::H || a.element == El::D;
70 atom.
aniso = {0, 0, 0, 0, 0, 0};
81 atom.
aniso = {u, u, u, 0.f, 0.f, 0.f};
100 for (
Residue& res : chain.residues) {
101 for (
Atom& atom : res.atoms)
102 atom.serial = ++serial;
119template<
typename Func>
145template<
typename Func>
148 for (
ModRes& modres :
st.mod_residues)
149 func(modres.chain_name, modres.res_id.seqid);
153 func(
sel.chain,
sel.res_begin);
154 func(
sel.chain,
sel.res_end);
160 auto update = [&](std::string& name) {
165 for (
ModRes& modres :
st.mod_residues)
166 update(modres.chain_name);
172 for (
Chain& chain : model.chains)
183 for (
ModRes& modres :
st.mod_residues)
186 for (std::string&
mon_ids :
ent.full_sequence)
187 for (
size_t start = 0;;) {
188 size_t end =
mon_ids.find(
',', start);
191 if (end != std::string::npos)
194 if (end == std::string::npos)
199 for (
Chain& chain : model.chains)
200 for (
Residue& res : chain.residues)
206 const std::map<std::string, std::string>&
old_new) {
217 for (
Chain& chain : model.chains)
218 for (
Residue& res : chain.residues)
220 for (
Atom& atom : res.atoms)
227 for (
size_t i = res.
atoms.size();
i-- != 0; ) {
233 if (atom.
name[0] ==
'H')
250 if (
deut->name[0] ==
'H')
259 for (
auto d = res.
atoms.end();
d-- != res.
atoms.begin(); )
260 if (
d->element ==
El::D) {
262 auto h = res.
atoms.begin();
263 for (;
h != res.
atoms.end(); ++
h)
264 if (
h->element ==
El::H &&
h->pos.approx(
d->pos, 1
e-9))
266 if (
h != res.
atoms.end()) {
268 h->fraction =
h->occ > 0.f ?
d->occ /
h->occ : 0.f;
272 if (
i !=
d &&
i !=
h && (
i->name ==
h->name ||
i->name ==
d->name))
296 st.has_d_fraction =
false;
298 for (
Chain& chain : model.chains)
299 for (
Residue& res : chain.residues)
302 st.has_d_fraction =
true;
310 if (!
st.cell.explicit_matrices || !
st.cell.is_crystal())
313 st.cell.explicit_matrices =
false;
314 st.cell.calculate_properties();
Data structures to store macromolecular structure models.
void assign_serial_numbers(Model &model, bool numbered_ter=false)
set atom site serial numbers to 1, 2, ..., optionally leaving gaps for TERs
void remove_anisou(T &obj)
Remove anisotropic ADP.
void remove_alternative_conformations(T &obj)
Remove alternative conformations.
void replace_d_fraction_with_altlocs(Residue &res)
void vector_remove_if(std::vector< T > &v, F &&condition)
bool replace_deuterium_with_fraction(Residue &res)
constexpr T clamp(T v, T lo, T hi)
void process_addresses(Structure &st, Func func)
Helper function for processing (usually: changing) names and numbers in AtomAddress instances in meta...
void rename_chain(Structure &st, const std::string &old_name, const std::string &new_name)
void rename_residues(Structure &st, const std::string &old_name, const std::string &new_name)
void transform_pos_and_adp(T &obj, const Transform &tr)
apply Transform to both atom's position and ADP
void remove_hydrogens(T &obj)
Remove hydrogens.
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,...
void rename_atom_names(Structure &st, const std::string &res_name, const std::map< std::string, std::string > &old_new)
void standardize_crystal_frame(Structure &st)
Convert coordinates to the standard coordinate system for the unit cell.
void assign_b_iso(T &obj, float b_min, float b_max)
Set isotropic ADP to the range (b_min, b_max).
void process_sequence_ids(Structure &st, Func func)
Takes func(const std::string& chain_name, gemmi::SeqId& seqid).
void ensure_anisou(T &obj)
Set absent ANISOU to value from B_iso.
constexpr double u_to_b()
Represents atom site in macromolecular structure (~100 bytes).
std::vector< Residue > residues
std::vector< Chain > chains
Non-crystallographic symmetry operation (such as in the MTRIXn record)
Coordinates in Angstroms - orthogonal (Cartesian) coordinates.
std::vector< Atom > atoms
SMat33< Real > transformed_by(const Mat33 &m) const
Utilities. Mostly for working with strings and vectors.