Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2013-2023 The plumed team
3 : (see the PEOPLE file at the root of the distribution for a list of names)
4 :
5 : See http://www.plumed.org for more information.
6 :
7 : This file is part of plumed, version 2.
8 :
9 : plumed is free software: you can redistribute it and/or modify
10 : it under the terms of the GNU Lesser General Public License as published by
11 : the Free Software Foundation, either version 3 of the License, or
12 : (at your option) any later version.
13 :
14 : plumed is distributed in the hope that it will be useful,
15 : but WITHOUT ANY WARRANTY; without even the implied warranty of
16 : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 : GNU Lesser General Public License for more details.
18 :
19 : You should have received a copy of the GNU Lesser General Public License
20 : along with plumed. If not, see <http://www.gnu.org/licenses/>.
21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
22 :
23 : #define cutOffNB 0.60 // buffer distance for neighbour-lists
24 : #define cutOffDist 0.50 // cut off distance for non-bonded pairwise forces
25 : #define cutOnDist 0.32 // cut off distance for non-bonded pairwise forces
26 : #define cutOffNB2 cutOffNB*cutOffNB // squared buffer distance for neighbour-lists
27 : #define cutOffDist2 cutOffDist*cutOffDist
28 : #define cutOnDist2 cutOnDist*cutOnDist
29 : #define invswitch 1.0/((cutOffDist2-cutOnDist2)*(cutOffDist2-cutOnDist2)*(cutOffDist2-cutOnDist2))
30 : #define cutOffDist4 cutOffDist2*cutOffDist2
31 : #define cutMixed cutOffDist2*cutOffDist2*cutOffDist2 -3.*cutOffDist2*cutOffDist2*cutOnDist2
32 :
33 : #include <string>
34 : #include <fstream>
35 : #include <iterator>
36 : #include <sstream>
37 :
38 : #include "MetainferenceBase.h"
39 : #include "core/ActionRegister.h"
40 : #include "tools/Pbc.h"
41 : #include "tools/PDB.h"
42 : #include "tools/Torsion.h"
43 :
44 : namespace PLMD {
45 : namespace isdb {
46 :
47 : //+PLUMEDOC ISDB_COLVAR CS2BACKBONE
48 : /*
49 : Calculates the backbone chemical shifts for a protein.
50 :
51 : The functional form is that of CamShift \cite Kohlhoff:2009us. The chemical shift
52 : of the selected nuclei can be saved as components. Alternatively one can calculate either
53 : the CAMSHIFT score (useful as a collective variable \cite Granata:2013dk or as a scoring
54 : function \cite Robustelli:2010dn) or a \ref METAINFERENCE score (using DOSCORE).
55 : For these two latter cases experimental chemical shifts must be provided.
56 :
57 : CS2BACKBONE calculation can be relatively heavy because it often uses a large number of atoms, it can
58 : be run in parallel using MPI and \ref Openmp.
59 :
60 : As a general rule, when using \ref CS2BACKBONE or other experimental restraints it may be better to
61 : increase the accuracy of the constraint algorithm due to the increased strain on the bonded structure.
62 : In the case of GROMACS it is safer to use lincs-iter=2 and lincs-order=6.
63 :
64 : In general the system for which chemical shifts are calculated must be completely included in
65 : ATOMS and a TEMPLATE pdb file for the same atoms should be provided as well in the folder DATADIR.
66 : The system is made automatically whole unless NOPBC is used, in particular if the system is made
67 : by multiple chains it is usually better to use NOPBC and make the molecule whole \ref WHOLEMOLECULES
68 : selecting an appropriate order of the atoms. The pdb file is needed to the generate a simple topology of the protein.
69 : For histidine residues in protonation states different from D the HIE/HSE HIP/HSP name should be used.
70 : GLH and ASH can be used for the alternative protonation of GLU and ASP. Non-standard amino acids and other
71 : molecules are not yet supported, but in principle they can be named UNK. If multiple chains are present
72 : the chain identifier must be in the standard PDB format, together with the TER keyword at the end of each chain.
73 : Termini groups like ACE or NME should be removed from the TEMPLATE pdb because they are not recognized by
74 : CS2BACKBONE.
75 :
76 : In addition to a pdb file one needs to provide a list of chemical shifts to be calculated using one
77 : file per nucleus type (CAshifts.dat, CBshifts.dat, Cshifts.dat, Hshifts.dat, HAshifts.dat, Nshifts.dat),
78 : add only the files for the nuclei you need, but each file should include all protein residues.
79 : A chemical shift for a nucleus is calculated if a value greater than 0 is provided.
80 : For practical purposes the value can correspond to the experimental value.
81 : Residues numbers should match that used in the pdb file, but must be positive, so double check the pdb.
82 : The first and last residue of each chain should be preceded by a # character.
83 :
84 : \verbatim
85 : CAshifts.dat:
86 : #1 0.0
87 : 2 55.5
88 : 3 58.4
89 : .
90 : .
91 : #last 0.0
92 : #first of second chain
93 : .
94 : #last of second chain
95 : \endverbatim
96 :
97 : The default behavior is to store the values for the active nuclei in components (ca-#, cb-#,
98 : co-#, ha-#, hn-#, nh-# and expca-#, expcb-#, expco-#, expha-#, exphn-#, exp-nh#) with NOEXP it is possible
99 : to only store the back-calculated values, where # includes a chain and residue number.
100 :
101 : One additional file is always needed in the folder DATADIR: camshift.db. This file includes all the parameters needed to
102 : calculate the chemical shifts and can be found in regtest/isdb/rt-cs2backbone/data/ .
103 :
104 : Additional material and examples can be also found in the tutorial \ref isdb-1 as well as in the cs2backbone regtests
105 : in the isdb folder.
106 :
107 : \par Examples
108 :
109 : In this first example the chemical shifts are used to calculate a collective variable to be used
110 : in NMR driven Metadynamics \cite Granata:2013dk :
111 :
112 : \plumedfile
113 : #SETTINGS AUXFOLDER=regtest/isdb/rt-cs2backbone/data
114 : whole: GROUP ATOMS=2612-2514:-1,961-1:-1,2466-962:-1,2513-2467:-1
115 : WHOLEMOLECULES ENTITY0=whole
116 : cs: CS2BACKBONE ATOMS=1-2612 DATADIR=data/ TEMPLATE=template.pdb CAMSHIFT NOPBC
117 : metad: METAD ARG=cs HEIGHT=0.5 SIGMA=0.1 PACE=200 BIASFACTOR=10
118 : PRINT ARG=cs,metad.bias FILE=COLVAR STRIDE=100
119 : \endplumedfile
120 :
121 : In this second example the chemical shifts are used as replica-averaged restrained as in \cite Camilloni:2012je \cite Camilloni:2013hs .
122 :
123 : \plumedfile
124 : #SETTINGS AUXFOLDER=regtest/isdb/rt-cs2backbone/data NREPLICAS=2
125 : cs: CS2BACKBONE ATOMS=1-174 DATADIR=data/
126 : encs: ENSEMBLE ARG=(cs\.hn-.*),(cs\.nh-.*)
127 : stcs: STATS ARG=encs.* SQDEVSUM PARARG=(cs\.exphn-.*),(cs\.expnh-.*)
128 : RESTRAINT ARG=stcs.sqdevsum AT=0 KAPPA=0 SLOPE=24
129 :
130 : PRINT ARG=(cs\.hn-.*),(cs\.nh-.*) FILE=RESTRAINT STRIDE=100
131 :
132 : \endplumedfile
133 :
134 : This third example show how to use chemical shifts to calculate a \ref METAINFERENCE score .
135 :
136 : \plumedfile
137 : #SETTINGS AUXFOLDER=regtest/isdb/rt-cs2backbone/data
138 : cs: CS2BACKBONE ATOMS=1-174 DATADIR=data/ SIGMA_MEAN0=1.0 DOSCORE
139 : csbias: BIASVALUE ARG=cs.score
140 :
141 : PRINT ARG=(cs\.hn-.*),(cs\.nh-.*) FILE=CS.dat STRIDE=1000
142 : PRINT ARG=cs.score FILE=BIAS STRIDE=100
143 : \endplumedfile
144 :
145 : */
146 : //+ENDPLUMEDOC
147 :
148 : class CS2BackboneDB {
149 : enum { STD, GLY, PRO};
150 : enum { HA_ATOM, H_ATOM, N_ATOM, CA_ATOM, CB_ATOM, C_ATOM };
151 : static const unsigned aa_kind = 3;
152 : static const unsigned atm_kind = 6;
153 : static const unsigned numXtraDists = 27;
154 :
155 : // ALA, ARG, ASN, ASP, CYS, GLU, GLN, GLY, HIS, ILE, LEU, LYS, MET, PHE, PRO, SER, THR, TRP, TYR, VAL
156 : double c_aa[aa_kind][atm_kind][20];
157 : double c_aa_prev[aa_kind][atm_kind][20];
158 : double c_aa_succ[aa_kind][atm_kind][20];
159 : double co_bb[aa_kind][atm_kind][16];
160 : double co_sc_[aa_kind][atm_kind][20][20];
161 : double co_xd[aa_kind][atm_kind][numXtraDists];
162 : double co_sphere[aa_kind][atm_kind][2][8];
163 : // for ring current effects
164 : // Phe, Tyr, Trp_1, Trp_2, His
165 : double co_ring[aa_kind][atm_kind][5];
166 : // for dihedral angles
167 : // co * (a * cos(3 * omega + c) + b * cos(omega + d))
168 : double co_da[aa_kind][atm_kind][3];
169 : double pars_da[aa_kind][atm_kind][3][5];
170 :
171 : public:
172 :
173 10926 : inline unsigned kind(const std::string &s) {
174 10926 : if(s=="GLY") return GLY;
175 9324 : else if(s=="PRO") return PRO;
176 : return STD;
177 : }
178 :
179 10926 : inline unsigned atom_kind(const std::string &s) {
180 10926 : if(s=="HA")return HA_ATOM;
181 10872 : else if(s=="H") return H_ATOM;
182 8010 : else if(s=="N") return N_ATOM;
183 5148 : else if(s=="CA")return CA_ATOM;
184 2160 : else if(s=="CB")return CB_ATOM;
185 54 : else if(s=="C") return C_ATOM;
186 : return -1;
187 : }
188 :
189 : unsigned get_numXtraDists() {return numXtraDists;}
190 :
191 : //PARAMETERS
192 : inline double * CONSTAACURR(const unsigned a_kind, const unsigned at_kind) {return c_aa[a_kind][at_kind];}
193 : inline double * CONSTAANEXT(const unsigned a_kind, const unsigned at_kind) {return c_aa_succ[a_kind][at_kind];}
194 : inline double * CONSTAAPREV(const unsigned a_kind, const unsigned at_kind) {return c_aa_prev[a_kind][at_kind];}
195 : inline double * CONST_BB2(const unsigned a_kind, const unsigned at_kind) {return co_bb[a_kind][at_kind];}
196 : inline double * CONST_SC2(const unsigned a_kind, const unsigned at_kind, unsigned res_type) { return co_sc_[a_kind][at_kind][res_type];}
197 : inline double * CONST_XD(const unsigned a_kind, const unsigned at_kind) { return co_xd[a_kind][at_kind];}
198 : inline double * CO_SPHERE(const unsigned a_kind, const unsigned at_kind, unsigned exp_type) { return co_sphere[a_kind][at_kind][exp_type];}
199 : inline double * CO_RING(const unsigned a_kind, const unsigned at_kind) { return co_ring[a_kind][at_kind];}
200 : inline double * CO_DA(const unsigned a_kind, const unsigned at_kind) { return co_da[a_kind][at_kind];}
201 : inline double * PARS_DA(const unsigned a_kind, const unsigned at_kind, const unsigned ang_kind) { return pars_da[a_kind][at_kind][ang_kind];}
202 :
203 18 : void parse(const std::string &file, const double dscale) {
204 18 : std::ifstream in;
205 18 : in.open(file.c_str());
206 18 : if(!in) plumed_merror("Unable to open DB file: " + file);
207 :
208 : unsigned c_kind = 0;
209 : unsigned c_atom = 0;
210 : unsigned nline = 0;
211 :
212 396 : for(unsigned i=0; i<3; i++) for(unsigned j=0; j<6; j++) {
213 6804 : for(unsigned k=0; k<20; k++) {
214 6480 : c_aa[i][j][k]=0.;
215 6480 : c_aa_prev[i][j][k]=0.;
216 6480 : c_aa_succ[i][j][k]=0.;
217 136080 : for(unsigned m=0; m<20; m++) co_sc_[i][j][k][m]=0.;
218 : }
219 5508 : for(unsigned k=0; k<16; k++) {co_bb[i][j][k]=0.; }
220 2916 : for(unsigned k=0; k<8; k++) { co_sphere[i][j][0][k]=0.; co_sphere[i][j][1][k]=0.; }
221 1296 : for(unsigned k=0; k<3; k++) {
222 972 : co_da[i][j][k]=0.;
223 5832 : for(unsigned l=0; l<5; l++) pars_da[i][j][k][l]=0.;
224 : }
225 1944 : for(unsigned k=0; k<5; k++) co_ring[i][j][k]=0.;
226 9072 : for(unsigned k=0; k<numXtraDists; k++) co_xd[i][j][k]=0.;
227 : }
228 :
229 37710 : while(!in.eof()) {
230 : std::string line;
231 37692 : getline(in,line);
232 : ++nline;
233 37692 : if(line.compare(0,1,"#")==0) continue;
234 : std::vector<std::string> tok;
235 : std::vector<std::string> tmp;
236 40860 : tok = split(line,' ');
237 261828 : for(unsigned q=0; q<tok.size(); q++)
238 241398 : if(tok[q].size()) tmp.push_back(tok[q]);
239 20430 : tok = tmp;
240 20430 : if(tok.size()==0) continue;
241 20412 : if(tok[0]=="PAR") {
242 324 : c_kind = kind(tok[2]);
243 324 : c_atom = atom_kind(tok[1]);
244 324 : continue;
245 : }
246 20088 : else if(tok[0]=="WEIGHT") {
247 324 : continue;
248 : }
249 19764 : else if(tok[0]=="FLATBTM") {
250 324 : continue;
251 : }
252 19440 : else if (tok[0] == "SCALEHARM") {
253 324 : continue;
254 : }
255 19116 : else if (tok[0] == "TANHAMPLI") {
256 324 : continue;
257 : }
258 18792 : else if (tok[0] == "ENDHARMON") {
259 324 : continue;
260 : }
261 18468 : else if (tok[0] == "MAXRCDEVI") {
262 324 : continue;
263 : }
264 18144 : else if (tok[0] == "RANDCOIL") {
265 324 : continue;
266 : }
267 17820 : else if (tok[0] == "CONST") {
268 324 : continue;
269 : }
270 17496 : else if (tok[0] == "CONSTAA") {
271 324 : assign(c_aa[c_kind][c_atom],tok,1);
272 324 : continue;
273 : }
274 17172 : else if (tok[0] == "CONSTAA-1") {
275 324 : assign(c_aa_prev[c_kind][c_atom],tok,1);
276 324 : continue;
277 : }
278 16848 : else if (tok[0] == "CONSTAA+1") {
279 324 : assign(c_aa_succ[c_kind][c_atom],tok,1);
280 324 : continue;
281 : }
282 16524 : else if (tok[0] == "COBB1") {
283 324 : continue;
284 : }
285 16200 : else if (tok[0] == "COBB2") {
286 : //angstrom to nm
287 324 : assign(co_bb[c_kind][c_atom],tok,dscale);
288 324 : continue;
289 : }
290 15876 : else if (tok[0] == "SPHERE1") {
291 : // angstrom^-3 to nm^-3
292 324 : assign(co_sphere[c_kind][c_atom][0],tok,1./(dscale*dscale*dscale));
293 324 : continue;
294 : }
295 15552 : else if (tok[0] == "SPHERE2") {
296 : //angstrom to nm
297 324 : assign(co_sphere[c_kind][c_atom][1],tok,dscale);
298 324 : continue;
299 : }
300 15228 : else if (tok[0] == "DIHEDRALS") {
301 324 : assign(co_da[c_kind][c_atom],tok,1);
302 324 : continue;
303 : }
304 14904 : else if (tok[0] == "RINGS") {
305 : // angstrom^-3 to nm^-3
306 324 : assign(co_ring[c_kind][c_atom],tok,1./(dscale*dscale*dscale));
307 1944 : for(unsigned i=1; i<tok.size(); i++)
308 1620 : co_ring[c_kind][c_atom][i-1] *= 1000;
309 324 : continue;
310 324 : }
311 14580 : else if (tok[0] == "HBONDS") {
312 324 : continue;
313 : }
314 14256 : else if (tok[0] == "XTRADISTS") {
315 : //angstrom to nm
316 324 : assign(co_xd[c_kind][c_atom],tok,dscale);
317 324 : continue;
318 : }
319 13932 : else if(tok[0]=="DIHEDPHI") {
320 324 : assign(pars_da[c_kind][c_atom][0],tok,1);
321 324 : continue;
322 : }
323 13608 : else if(tok[0]=="DIHEDPSI") {
324 324 : assign(pars_da[c_kind][c_atom][1],tok,1);
325 324 : continue;
326 : }
327 13284 : else if(tok[0]=="DIHEDCHI1") {
328 324 : assign(pars_da[c_kind][c_atom][2],tok,1);
329 324 : continue;
330 : }
331 :
332 : bool ok = false;
333 : std::string scIdent1 [] = {"COSCALA1", "COSCARG1", "COSCASN1", "COSCASP1", "COSCCYS1", "COSCGLN1", "COSCGLU1",
334 : "COSCGLY1", "COSCHIS1", "COSCILE1", "COSCLEU1", "COSCLYS1", "COSCMET1", "COSCPHE1",
335 : "COSCPRO1", "COSCSER1", "COSCTHR1", "COSCTRP1", "COSCTYR1", "COSCVAL1"
336 272160 : };
337 :
338 204120 : for(unsigned scC = 0; scC < 20; scC++) {
339 197640 : if(tok[0]==scIdent1[scC]) {
340 : ok = true;
341 : break;
342 : }
343 : }
344 285120 : if(ok) continue;
345 :
346 : std::string scIdent2 [] = {"COSCALA2", "COSCARG2", "COSCASN2", "COSCASP2", "COSCCYS2", "COSCGLN2", "COSCGLU2",
347 : "COSCGLY2", "COSCHIS2", "COSCILE2", "COSCLEU2", "COSCLYS2", "COSCMET2", "COSCPHE2",
348 : "COSCPRO2", "COSCSER2", "COSCTHR2", "COSCTRP2", "COSCTYR2", "COSCVAL2"
349 136080 : };
350 :
351 68040 : for(unsigned scC = 0; scC < 20; scC++) {
352 68040 : if(tok[0]==scIdent2[scC]) {
353 : //angstrom to nm
354 6480 : assign(co_sc_[c_kind][c_atom][scC],tok,dscale);
355 : ok = true; break;
356 : }
357 : }
358 142560 : if(ok) continue;
359 :
360 0 : if(tok.size()) {
361 0 : std::string str_err = "DB WARNING: unrecognized token: " + tok[0];
362 0 : plumed_merror(str_err);
363 : }
364 428670 : }
365 18 : in.close();
366 18 : }
367 :
368 : private:
369 :
370 20430 : std::vector<std::string> &split(const std::string &s, char delim, std::vector<std::string> &elems) {
371 20430 : std::stringstream ss(s);
372 : std::string item;
373 261828 : while (getline(ss, item, delim)) {
374 241398 : elems.push_back(item);
375 : }
376 20430 : return elems;
377 20430 : }
378 :
379 20430 : std::vector<std::string> split(const std::string &s, char delim) {
380 : std::vector<std::string> elems;
381 20430 : split(s, delim, elems);
382 20430 : return elems;
383 0 : }
384 :
385 10368 : void assign(double * f, const std::vector<std::string> & v, const double scale) {
386 122796 : for(unsigned i=1; i<v.size(); i++) {
387 112428 : f[i-1] = scale*(atof(v[i].c_str()));
388 112428 : if(std::abs(f[i-1])<0.000001) f[i-1]=0.;
389 : }
390 10368 : }
391 : };
392 :
393 : class CS2Backbone : public MetainferenceBase {
394 : struct ChemicalShift {
395 : double exp_cs; // a reference chemical shifts
396 : Value *comp; // a pointer to the component
397 : unsigned res_kind; // residue type (STD/GLY/PRO)
398 : unsigned atm_kind; // nuclues (HA/CA/CB/CO/NH/HN)
399 : unsigned res_type_prev; // previous residue (ALA/VAL/..)
400 : unsigned res_type_curr; // current residue (ALA/VAL/..)
401 : unsigned res_type_next; // next residue (ALA/VAL/..)
402 : std::string res_name; // residue name
403 : std::string nucleus; // chemical shift
404 : bool has_chi1; // does we have a chi1
405 : unsigned csatoms; // fixed number of atoms used
406 : unsigned totcsatoms; // number of atoms used
407 : unsigned res_num; // residue number
408 : unsigned chain; // chain number
409 : unsigned ipos; // index of the atom for which we are calculating the chemical shifts
410 : std::vector<unsigned> bb; // atoms for the previous, current and next backbone
411 : std::vector<unsigned> side_chain;// atoms for the current sidechain
412 : std::vector<int> xd1; // additional couple of atoms
413 : std::vector<int> xd2; // additional couple of atoms
414 : std::vector<unsigned> box_nb; // non-bonded atoms
415 :
416 10602 : ChemicalShift():
417 10602 : exp_cs(0.),
418 10602 : comp(NULL),
419 10602 : res_kind(0),
420 10602 : atm_kind(0),
421 10602 : res_type_prev(0),
422 10602 : res_type_curr(0),
423 10602 : res_type_next(0),
424 10602 : res_name(""),
425 10602 : nucleus(""),
426 10602 : has_chi1(true),
427 10602 : csatoms(0),
428 10602 : totcsatoms(0),
429 10602 : res_num(0),
430 10602 : chain(0),
431 10602 : ipos(0)
432 : {
433 10602 : xd1.reserve(26);
434 10602 : xd2.reserve(26);
435 10602 : box_nb.reserve(150);
436 10602 : }
437 : };
438 :
439 : struct RingInfo {
440 : enum {R_PHE, R_TYR, R_TRP1, R_TRP2, R_HIS};
441 : unsigned rtype; // one out of five different types
442 : unsigned atom[6]; // up to six member per ring
443 : unsigned numAtoms; // number of ring members (5 or 6)
444 : Vector position; // center of ring coordinates
445 : Vector normVect; // ring plane normal std::vector
446 : Vector g[6]; // std::vector of the std::vectors used for normVect
447 : double lengthN2; // square of length of normVect
448 : double lengthNV; // length of normVect
449 360 : RingInfo():
450 360 : rtype(0),
451 360 : numAtoms(0),
452 360 : lengthN2(NAN),
453 2520 : lengthNV(NAN)
454 : {
455 2520 : for(unsigned i=0; i<6; i++) atom[i]=0;
456 360 : }
457 : };
458 :
459 : enum aa_t {ALA, ARG, ASN, ASP, CYS, GLN, GLU, GLY, HIS, ILE, LEU, LYS, MET, PHE, PRO, SER, THR, TRP, TYR, VAL, UNK};
460 : enum sequence_t {Np, CAp, HAp, Cp, Op, Nc, Hc, CAc, HAc, Cc, Oc, Nn, Hn, CAn, HAn, Cn, CBc, CGc};
461 :
462 : CS2BackboneDB db;
463 : std::vector<ChemicalShift> chemicalshifts;
464 :
465 : std::vector<RingInfo> ringInfo;
466 : std::vector<unsigned> type;
467 : std::vector<unsigned> res_num;
468 : unsigned max_cs_atoms;
469 : unsigned box_nupdate;
470 : unsigned box_count;
471 : bool camshift;
472 : bool pbc;
473 : bool serial;
474 :
475 : void init_cs(const std::string &file, const std::string &k, const PDB &pdb);
476 : void update_neighb();
477 : void compute_ring_parameters();
478 : void init_types(const PDB &pdb);
479 : void init_rings(const PDB &pdb);
480 : aa_t frag2enum(const std::string &aa);
481 : std::vector<std::string> side_chain_atoms(const std::string &s);
482 : bool isSP2(const std::string & resType, const std::string & atomName);
483 : bool is_chi1_cx(const std::string & frg, const std::string & atm);
484 : void xdist_name_map(std::string & name);
485 :
486 : public:
487 :
488 : explicit CS2Backbone(const ActionOptions&);
489 : static void registerKeywords( Keywords& keys );
490 : void calculate() override;
491 : void update() override;
492 : };
493 :
494 10455 : PLUMED_REGISTER_ACTION(CS2Backbone,"CS2BACKBONE")
495 :
496 19 : void CS2Backbone::registerKeywords( Keywords& keys ) {
497 19 : componentsAreNotOptional(keys);
498 19 : MetainferenceBase::registerKeywords( keys );
499 38 : keys.addFlag("NOPBC",false,"ignore the periodic boundary conditions when calculating distances");
500 38 : keys.addFlag("SERIAL",false,"Perform the calculation in serial - for debug purpose");
501 38 : keys.add("atoms","ATOMS","The atoms to be included in the calculation, e.g. the whole protein.");
502 38 : keys.add("compulsory","DATADIR","data/","The folder with the experimental chemical shifts.");
503 38 : keys.add("compulsory","TEMPLATE","template.pdb","A PDB file of the protein system.");
504 38 : keys.add("compulsory","NEIGH_FREQ","20","Period in step for neighbor list update.");
505 38 : keys.addFlag("CAMSHIFT",false,"Set to TRUE if you to calculate a single CamShift score.");
506 38 : keys.addFlag("NOEXP",false,"Set to TRUE if you don't want to have fixed components with the experimental values.");
507 38 : keys.addOutputComponent("ha","default","the calculated Ha hydrogen chemical shifts");
508 38 : keys.addOutputComponent("hn","default","the calculated H hydrogen chemical shifts");
509 38 : keys.addOutputComponent("nh","default","the calculated N nitrogen chemical shifts");
510 38 : keys.addOutputComponent("ca","default","the calculated Ca carbon chemical shifts");
511 38 : keys.addOutputComponent("cb","default","the calculated Cb carbon chemical shifts");
512 38 : keys.addOutputComponent("co","default","the calculated C' carbon chemical shifts");
513 38 : keys.addOutputComponent("expha","default","the experimental Ha hydrogen chemical shifts");
514 38 : keys.addOutputComponent("exphn","default","the experimental H hydrogen chemical shifts");
515 38 : keys.addOutputComponent("expnh","default","the experimental N nitrogen chemical shifts");
516 38 : keys.addOutputComponent("expca","default","the experimental Ca carbon chemical shifts");
517 38 : keys.addOutputComponent("expcb","default","the experimental Cb carbon chemical shifts");
518 38 : keys.addOutputComponent("expco","default","the experimental C' carbon chemical shifts");
519 19 : }
520 :
521 18 : CS2Backbone::CS2Backbone(const ActionOptions&ao):
522 : PLUMED_METAINF_INIT(ao),
523 18 : max_cs_atoms(0),
524 18 : camshift(false),
525 18 : pbc(true),
526 18 : serial(false)
527 : {
528 : std::vector<AtomNumber> used_atoms;
529 18 : parseAtomList("ATOMS",used_atoms);
530 :
531 18 : parseFlag("CAMSHIFT",camshift);
532 18 : if(camshift&&getDoScore()) plumed_merror("It is not possible to use CAMSHIFT and DOSCORE at the same time");
533 :
534 18 : bool nopbc=!pbc;
535 18 : parseFlag("NOPBC",nopbc);
536 18 : pbc=!nopbc;
537 :
538 18 : parseFlag("SERIAL",serial);
539 :
540 18 : bool noexp=false;
541 36 : parseFlag("NOEXP",noexp);
542 :
543 : std::string stringa_data;
544 36 : parse("DATADIR",stringa_data);
545 :
546 : std::string stringa_template;
547 18 : parse("TEMPLATE",stringa_template);
548 :
549 18 : box_count=0;
550 18 : box_nupdate=20;
551 18 : parse("NEIGH_FREQ", box_nupdate);
552 :
553 18 : std::string stringadb = stringa_data + std::string("/camshift.db");
554 36 : std::string stringapdb = stringa_data + std::string("/") + stringa_template;
555 :
556 : /* Length conversion (parameters are tuned for angstrom) */
557 : double scale=1.;
558 18 : if(!plumed.getAtoms().usingNaturalUnits()) {
559 18 : scale = 10.*atoms.getUnits().getLength();
560 : }
561 :
562 18 : log.printf(" Initialization of the predictor ...\n");
563 18 : db.parse(stringadb,scale);
564 :
565 18 : PDB pdb;
566 36 : if( !pdb.read(stringapdb,plumed.getAtoms().usingNaturalUnits(),1./scale) ) plumed_merror("missing input file " + stringapdb);
567 :
568 : // first of all we build the list of chemical shifts we want to predict
569 18 : log.printf(" Reading experimental data ...\n"); log.flush();
570 36 : stringadb = stringa_data + std::string("/CAshifts.dat");
571 18 : log.printf(" Initializing CA shifts %s\n", stringadb.c_str());
572 18 : init_cs(stringadb, "CA", pdb);
573 36 : stringadb = stringa_data + std::string("/CBshifts.dat");
574 18 : log.printf(" Initializing CB shifts %s\n", stringadb.c_str());
575 18 : init_cs(stringadb, "CB", pdb);
576 36 : stringadb = stringa_data + std::string("/Cshifts.dat");
577 18 : log.printf(" Initializing C' shifts %s\n", stringadb.c_str());
578 18 : init_cs(stringadb, "C", pdb);
579 36 : stringadb = stringa_data + std::string("/HAshifts.dat");
580 18 : log.printf(" Initializing HA shifts %s\n", stringadb.c_str());
581 18 : init_cs(stringadb, "HA", pdb);
582 36 : stringadb = stringa_data + std::string("/Hshifts.dat");
583 18 : log.printf(" Initializing H shifts %s\n", stringadb.c_str());
584 18 : init_cs(stringadb, "H", pdb);
585 36 : stringadb = stringa_data + std::string("/Nshifts.dat");
586 18 : log.printf(" Initializing N shifts %s\n", stringadb.c_str());
587 36 : init_cs(stringadb, "N", pdb);
588 :
589 18 : if(chemicalshifts.size()==0) plumed_merror("There are no chemical shifts to calculate, there must be at least a not empty file (CA|CB|C|HA|H|N|shifts.dat)");
590 :
591 18 : init_types(pdb);
592 18 : init_rings(pdb);
593 :
594 18 : log<<" Bibliography "
595 36 : <<plumed.cite("Kohlhoff K, Robustelli P, Cavalli A, Salvatella A, Vendruscolo M, J. Am. Chem. Soc. 131, 13894 (2009)");
596 28 : if(camshift) log<<plumed.cite("Granata D, Camilloni C, Vendruscolo M, Laio A, Proc. Natl. Acad. Sci. USA 110, 6817 (2013)");
597 39 : else log<<plumed.cite("Camilloni C, Robustelli P, De Simone A, Cavalli A, Vendruscolo M, J. Am. Chem. Soc. 134, 3968 (2012)");
598 36 : log<<plumed.cite("Bonomi M, Camilloni C, Bioinformatics, 33, 3999 (2017)");
599 18 : log<<"\n";
600 :
601 18 : if(camshift) {
602 5 : noexp = true;
603 5 : addValueWithDerivatives();
604 5 : setNotPeriodic();
605 : } else {
606 7670 : for(unsigned cs=0; cs<chemicalshifts.size(); cs++) {
607 7657 : std::string num; Tools::convert(chemicalshifts[cs].res_num,num);
608 7657 : std::string chain_num; Tools::convert(chemicalshifts[cs].chain,chain_num);
609 7657 : if(getDoScore()) {
610 4712 : addComponent(chemicalshifts[cs].nucleus+chain_num+"-"+num);
611 4712 : componentIsNotPeriodic(chemicalshifts[cs].nucleus+chain_num+"-"+num);
612 2356 : chemicalshifts[cs].comp = getPntrToComponent(chemicalshifts[cs].nucleus+chain_num+"-"+num);
613 4712 : setParameter(chemicalshifts[cs].exp_cs);
614 : } else {
615 10602 : addComponentWithDerivatives(chemicalshifts[cs].nucleus+chain_num+"-"+num);
616 10602 : componentIsNotPeriodic(chemicalshifts[cs].nucleus+chain_num+"-"+num);
617 5301 : chemicalshifts[cs].comp = getPntrToComponent(chemicalshifts[cs].nucleus+chain_num+"-"+num);
618 : }
619 : }
620 13 : if(getDoScore()) Initialise(chemicalshifts.size());
621 : }
622 :
623 18 : if(!noexp) {
624 7670 : for(unsigned cs=0; cs<chemicalshifts.size(); cs++) {
625 7657 : std::string num; Tools::convert(chemicalshifts[cs].res_num,num);
626 7657 : std::string chain_num; Tools::convert(chemicalshifts[cs].chain,chain_num);
627 15314 : addComponent("exp"+chemicalshifts[cs].nucleus+chain_num+"-"+num);
628 15314 : componentIsNotPeriodic("exp"+chemicalshifts[cs].nucleus+chain_num+"-"+num);
629 15314 : Value* comp=getPntrToComponent("exp"+chemicalshifts[cs].nucleus+chain_num+"-"+num);
630 7657 : comp->set(chemicalshifts[cs].exp_cs);
631 : }
632 : }
633 :
634 18 : requestAtoms(used_atoms, false);
635 18 : setDerivatives();
636 18 : checkRead();
637 36 : }
638 :
639 108 : void CS2Backbone::init_cs(const std::string &file, const std::string &nucl, const PDB &pdb) {
640 : // number of chains
641 : std::vector<std::string> chains;
642 108 : pdb.getChainNames( chains );
643 : unsigned ichain=0;
644 :
645 108 : std::ifstream in;
646 108 : in.open(file.c_str());
647 108 : if(!in) return;
648 108 : std::istream_iterator<std::string> iter(in), end;
649 : unsigned begin=0;
650 :
651 108 : while(iter!=end) {
652 19008 : std::string tok = *iter;
653 : ++iter;
654 19008 : if(tok[0]=='#') {
655 : ++iter;
656 432 : if(begin==1) {
657 : begin=0;
658 216 : ichain++;
659 : } else begin=1;
660 432 : continue;
661 : }
662 : int ro = atoi(tok.c_str());
663 18576 : if(ro<0) plumed_merror("Residue numbers should be positive\n");
664 18576 : unsigned resnum = static_cast<unsigned> (ro);
665 : tok = *iter;
666 : ++iter;
667 : double cs = atof(tok.c_str());
668 18576 : if(cs==0) continue;
669 :
670 : unsigned fres, lres;
671 : std::string errmsg;
672 11070 : pdb.getResidueRange(chains[ichain], fres, lres, errmsg);
673 11070 : if(resnum==fres||resnum==lres) plumed_merror("First and Last residue of each chain should be annotated as # in " + file + " Remember that residue numbers should match");
674 :
675 : // check in the PDB for the chain/residue/atom and enable the chemical shift
676 11070 : std::string RES = pdb.getResidueName(resnum, chains[ichain]);
677 98766 : if(RES=="HIE"||RES=="HIP"||RES=="HIS"||RES=="HSP"||RES=="HSE"||RES=="CYS"||RES=="GLH"||RES=="ASH"||RES=="UNK") continue;
678 10944 : if(RES=="GLN"&&nucl=="CB") continue;
679 11322 : if(RES=="ILE"&&nucl=="CB") continue;
680 10998 : if(RES=="PRO"&&nucl=="N") continue;
681 10998 : if(RES=="PRO"&&nucl=="H") continue;
682 10998 : if(RES=="PRO"&&nucl=="CB") continue;
683 12240 : if(RES=="GLY"&&nucl=="HA") continue;
684 12240 : if(RES=="GLY"&&nucl=="CB") continue;
685 :
686 10602 : ChemicalShift tmp_cs;
687 :
688 10602 : tmp_cs.exp_cs = cs;
689 10602 : if(nucl=="CA") tmp_cs.nucleus = "ca-";
690 7668 : else if(nucl=="CB") tmp_cs.nucleus = "cb-";
691 5616 : else if(nucl=="C") tmp_cs.nucleus = "co-";
692 5616 : else if(nucl=="HA") tmp_cs.nucleus = "ha-";
693 5616 : else if(nucl=="H") tmp_cs.nucleus = "hn-";
694 2808 : else if(nucl=="N") tmp_cs.nucleus = "nh-";
695 10602 : tmp_cs.chain = ichain;
696 10602 : tmp_cs.res_num = resnum;
697 10602 : tmp_cs.res_type_curr = frag2enum(RES);
698 10602 : tmp_cs.res_type_prev = frag2enum(pdb.getResidueName(resnum-1, chains[ichain]));
699 10602 : tmp_cs.res_type_next = frag2enum(pdb.getResidueName(resnum+1, chains[ichain]));
700 : tmp_cs.res_name = RES;
701 10602 : tmp_cs.res_kind = db.kind(RES);
702 10602 : tmp_cs.atm_kind = db.atom_kind(nucl);
703 20556 : if(RES!="ALA"&&RES!="GLY") {tmp_cs.bb.resize(18); tmp_cs.has_chi1=true;}
704 2142 : else {tmp_cs.bb.resize(16); tmp_cs.has_chi1=false;}
705 :
706 10602 : std::vector<AtomNumber> res_atoms = pdb.getAtomsInResidue(resnum, chains[ichain]);
707 : // find the position of the nucleus and of the other backbone atoms as well as for phi/psi/chi
708 173268 : for(unsigned a=0; a<res_atoms.size(); a++) {
709 162666 : std::string AN = pdb.getAtomName(res_atoms[a]);
710 162666 : if(nucl=="HA"&&(AN=="HA"||AN=="HA1"||AN=="HA3")) tmp_cs.ipos = res_atoms[a].index();
711 244530 : else if(nucl=="H"&&(AN=="H"||AN=="HN")) tmp_cs.ipos = res_atoms[a].index();
712 202248 : else if(nucl=="N"&&AN=="N") tmp_cs.ipos = res_atoms[a].index();
713 200862 : else if(nucl=="CA"&&AN=="CA") tmp_cs.ipos = res_atoms[a].index();
714 188244 : else if(nucl=="CB"&&AN=="CB") tmp_cs.ipos = res_atoms[a].index();
715 152064 : else if(nucl=="C"&&AN=="C" ) tmp_cs.ipos = res_atoms[a].index();
716 : }
717 :
718 10602 : std::vector<AtomNumber> prev_res_atoms = pdb.getAtomsInResidue(resnum-1, chains[ichain]);
719 : // find the position of the previous residues backbone atoms
720 168498 : for(unsigned a=0; a<prev_res_atoms.size(); a++) {
721 157896 : std::string AN = pdb.getAtomName(prev_res_atoms[a]);
722 157896 : if(AN=="N") { tmp_cs.bb[Np] = prev_res_atoms[a].index(); }
723 147294 : else if(AN=="CA") { tmp_cs.bb[CAp] = prev_res_atoms[a].index(); }
724 390690 : else if(AN=="HA"||AN=="HA1"||AN=="HA3") { tmp_cs.bb[HAp] = prev_res_atoms[a].index(); }
725 126090 : else if(AN=="C" ) { tmp_cs.bb[Cp] = prev_res_atoms[a].index(); }
726 115488 : else if(AN=="O" ) { tmp_cs.bb[Op] = prev_res_atoms[a].index(); }
727 : }
728 :
729 173268 : for(unsigned a=0; a<res_atoms.size(); a++) {
730 162666 : std::string AN = pdb.getAtomName(res_atoms[a]);
731 162666 : if(AN=="N") { tmp_cs.bb[Nc] = res_atoms[a].index(); }
732 438282 : else if(AN=="H" ||AN=="HN"||(AN=="CD"&&RES=="PRO")) { tmp_cs.bb[Hc] = res_atoms[a].index(); }
733 141462 : else if(AN=="CA") { tmp_cs.bb[CAc] = res_atoms[a].index(); }
734 372870 : else if(AN=="HA"||AN=="HA1"||AN=="HA3") { tmp_cs.bb[HAc] = res_atoms[a].index(); }
735 120258 : else if(AN=="C" ) { tmp_cs.bb[Cc] = res_atoms[a].index(); }
736 109656 : else if(AN=="O" ) { tmp_cs.bb[Oc] = res_atoms[a].index(); }
737 :
738 318852 : if(RES!="ALA"&&RES!="GLY") {
739 145728 : if(AN=="CB") tmp_cs.bb[CBc] = res_atoms[a].index();
740 145728 : if(is_chi1_cx(RES,AN)) tmp_cs.bb[CGc] = res_atoms[a].index();
741 : }
742 : }
743 :
744 10602 : std::vector<AtomNumber> next_res_atoms = pdb.getAtomsInResidue(resnum+1, chains[ichain]);
745 10602 : std::string NRES = pdb.getResidueName(resnum+1, chains[ichain]);
746 : // find the position of the previous residues backbone atoms
747 168948 : for(unsigned a=0; a<next_res_atoms.size(); a++) {
748 158346 : std::string AN = pdb.getAtomName(next_res_atoms[a]);
749 158346 : if(AN=="N") { tmp_cs.bb[Nn] = next_res_atoms[a].index(); }
750 426096 : else if(AN=="H" ||AN=="HN"||(AN=="CD"&&NRES=="PRO")) { tmp_cs.bb[Hn] = next_res_atoms[a].index(); }
751 137142 : else if(AN=="CA") { tmp_cs.bb[CAn] = next_res_atoms[a].index(); }
752 360198 : else if(AN=="HA"||AN=="HA1"||AN=="HA3") { tmp_cs.bb[HAn] = next_res_atoms[a].index(); }
753 115938 : else if(AN=="C" ) { tmp_cs.bb[Cn] = next_res_atoms[a].index(); }
754 : }
755 :
756 : // set sidechain atoms
757 10602 : std::vector<std::string> sc_atm = side_chain_atoms(RES);
758 :
759 141084 : for(unsigned sc=0; sc<sc_atm.size(); sc++) {
760 2501208 : for(unsigned aa=0; aa<res_atoms.size(); aa++) {
761 2370726 : if(pdb.getAtomName(res_atoms[aa])==sc_atm[sc]) {
762 99162 : tmp_cs.side_chain.push_back(res_atoms[aa].index());
763 : }
764 : }
765 : }
766 :
767 : // find atoms for extra distances
768 286254 : const std::string atomsP1[] = {"H", "H", "H", "C", "C", "C", "O", "O", "O", "N", "N", "N", "O", "O", "O", "N", "N", "N", "CG", "CG", "CG", "CG", "CG", "CG", "CG", "CA"};
769 10602 : const int resOffsetP1[] = { 0, 0, 0, -1, -1, -1, 0, 0, 0, 1, 1, 1, -1, -1, -1, 0, 0, 0, 0, 0, 0, 0, 0, -1, 1, -1};
770 :
771 286254 : const std::string atomsP2[] = {"HA", "C", "CB", "HA", "C", "CB", "HA", "N", "CB", "HA", "N", "CB", "HA", "N", "CB", "HA", "N", "CB", "HA", "N", "C", "C", "N", "CA", "CA", "CA"};
772 10602 : const int resOffsetP2[] = { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, -1, -1, -1, -1, -1, 0, 0, 0, -1, 1, 0, 0, 1};
773 :
774 286254 : for(unsigned q=0; q<db.get_numXtraDists()-1; q++) {
775 : std::vector<AtomNumber> at1;
776 275652 : if(resOffsetP1[q]== 0) at1 = res_atoms;
777 275652 : if(resOffsetP1[q]==-1) at1 = prev_res_atoms;
778 275652 : if(resOffsetP1[q]==+1) at1 = next_res_atoms;
779 :
780 : std::vector<AtomNumber> at2;
781 275652 : if(resOffsetP2[q]== 0) at2 = res_atoms;
782 275652 : if(resOffsetP2[q]==-1) at2 = prev_res_atoms;
783 275652 : if(resOffsetP2[q]==+1) at2 = next_res_atoms;
784 :
785 275652 : int tmp1 = -1;
786 2197566 : for(unsigned a=0; a<at1.size(); a++) {
787 2181708 : std::string name = pdb.getAtomName(at1[a]);
788 2181708 : xdist_name_map(name);
789 :
790 2181708 : if(name==atomsP1[q]) {
791 259794 : tmp1 = at1[a].index();
792 : break;
793 : }
794 : }
795 :
796 275652 : int tmp2 = -1;
797 1427688 : for(unsigned a=0; a<at2.size(); a++) {
798 1418076 : std::string name = pdb.getAtomName(at2[a]);
799 1418076 : xdist_name_map(name);
800 :
801 1418076 : if(name==atomsP2[q]) {
802 266040 : tmp2 = at2[a].index();
803 : break;
804 : }
805 : }
806 :
807 275652 : tmp_cs.xd1.push_back(tmp1);
808 275652 : tmp_cs.xd2.push_back(tmp2);
809 : }
810 :
811 : // ready to add a new chemical shifts
812 10602 : tmp_cs.csatoms = 1 + 16 + tmp_cs.side_chain.size() + 2*tmp_cs.xd1.size();
813 20556 : if(tmp_cs.res_name!="ALA"&&tmp_cs.res_name!="GLY") tmp_cs.csatoms += 2;
814 10602 : chemicalshifts.push_back(tmp_cs);
815 593712 : }
816 :
817 108 : in.close();
818 108 : }
819 :
820 : // this assigns an atom-type to each atom of the pdb
821 18 : void CS2Backbone::init_types(const PDB &pdb) {
822 : enum atom_t {D_C, D_H, D_N, D_O, D_S, D_C2, D_N2, D_O2};
823 18 : std::vector<AtomNumber> aa = pdb.getAtomNumbers();
824 47034 : for(unsigned i=0; i<aa.size(); i++) {
825 47016 : unsigned frag = pdb.getResidueNumber(aa[i]);
826 47016 : std::string fragName = pdb.getResidueName(aa[i]);
827 47016 : std::string atom_name = pdb.getAtomName(aa[i]);
828 47016 : char atom_type = atom_name[0];
829 47016 : if(isdigit(atom_name[0])) atom_type = atom_name[1];
830 47016 : res_num.push_back(frag);
831 47016 : unsigned t = 0;
832 47016 : if (!isSP2(fragName, atom_name)) {
833 : if (atom_type == 'C') t = D_C;
834 468 : else if (atom_type == 'O') t = D_O;
835 23256 : else if (atom_type == 'H') t = D_H;
836 4032 : else if (atom_type == 'N') t = D_N;
837 162 : else if (atom_type == 'S') t = D_S;
838 0 : else plumed_merror("Unknown atom type: " + atom_name);
839 : } else {
840 10080 : if (atom_type == 'C') t = D_C2;
841 4104 : else if (atom_type == 'O') t = D_O2;
842 0 : else if (atom_type == 'N') t = D_N2;
843 0 : else plumed_merror("Unknown atom type: " + atom_name);
844 : }
845 47016 : type.push_back(t);
846 : }
847 18 : }
848 :
849 18 : void CS2Backbone::init_rings(const PDB &pdb)
850 : {
851 126 : const std::string pheTyr_n[] = {"CG","CD1","CE1","CZ","CE2","CD2"};
852 126 : const std::string trp1_n[] = {"CD2","CE2","CZ2","CH2","CZ3","CE3"};
853 108 : const std::string trp2_n[] = {"CG","CD1","NE1","CE2","CD2"};
854 108 : const std::string his_n[] = {"CG","ND1","CD2","CE1","NE2"};
855 :
856 : // number of chains
857 : std::vector<std::string> chains;
858 18 : pdb.getChainNames( chains );
859 : unsigned total_rings_atoms = 0;
860 :
861 : // cycle over chains
862 54 : for(unsigned i=0; i<chains.size(); i++) {
863 : unsigned start, end;
864 : std::string errmsg;
865 36 : pdb.getResidueRange( chains[i], start, end, errmsg );
866 : // cycle over residues
867 3168 : for(unsigned res=start; res<end; res++) {
868 3132 : std::string frg = pdb.getResidueName(res, chains[i]);
869 14364 : if(!((frg=="PHE")||(frg=="TYR")||(frg=="TRP")||
870 8370 : (frg=="HIS")||(frg=="HIP")||(frg=="HID")||
871 5580 : (frg=="HIE")||(frg=="HSD")||(frg=="HSE")||
872 : (frg=="HSP"))) continue;
873 :
874 342 : std::vector<AtomNumber> frg_atoms = pdb.getAtomsInResidue(res,chains[i]);
875 :
876 396 : if(frg=="PHE"||frg=="TYR") {
877 324 : RingInfo ri;
878 6840 : for(unsigned a=0; a<frg_atoms.size(); a++) {
879 : unsigned atm = frg_atoms[a].index();
880 38808 : for(unsigned aa=0; aa<6; aa++) {
881 34236 : if(pdb.getAtomName(frg_atoms[a])==pheTyr_n[aa]) {
882 1944 : ri.atom[aa] = atm;
883 1944 : break;
884 : }
885 : }
886 : }
887 324 : ri.numAtoms = 6;
888 324 : total_rings_atoms += 6;
889 324 : if(frg=="PHE") ri.rtype = RingInfo::R_PHE;
890 324 : if(frg=="TYR") ri.rtype = RingInfo::R_TYR;
891 324 : ringInfo.push_back(ri);
892 :
893 18 : } else if(frg=="TRP") {
894 : //First ring
895 18 : RingInfo ri;
896 450 : for(unsigned a=0; a<frg_atoms.size(); a++) {
897 : unsigned atm = frg_atoms[a].index();
898 2646 : for(unsigned aa=0; aa<6; aa++) {
899 2322 : if(pdb.getAtomName(frg_atoms[a])==trp1_n[aa]) {
900 108 : ri.atom[aa] = atm;
901 108 : break;
902 : }
903 : }
904 : }
905 18 : ri.numAtoms = 6;
906 : total_rings_atoms += 6;
907 18 : ri.rtype = RingInfo::R_TRP1;
908 18 : ringInfo.push_back(ri);
909 : //Second Ring
910 18 : RingInfo ri2;
911 450 : for(unsigned a=0; a<frg_atoms.size(); a++) {
912 : unsigned atm = frg_atoms[a].index();
913 2322 : for(unsigned aa=0; aa<5; aa++) {
914 1980 : if(pdb.getAtomName(frg_atoms[a])==trp2_n[aa]) {
915 90 : ri2.atom[aa] = atm;
916 90 : break;
917 : }
918 : }
919 : }
920 18 : ri2.numAtoms = 5;
921 18 : total_rings_atoms += 3;
922 18 : ri2.rtype = RingInfo::R_TRP2;
923 18 : ringInfo.push_back(ri2);
924 :
925 0 : } else if((frg=="HIS")||(frg=="HIP")||(frg=="HID")||
926 0 : (frg=="HIE")||(frg=="HSD")||(frg=="HSE")||
927 : (frg=="HSP")) {//HIS case
928 0 : RingInfo ri;
929 0 : for(unsigned a=0; a<frg_atoms.size(); a++) {
930 : unsigned atm = frg_atoms[a].index();
931 0 : for(unsigned aa=0; aa<5; aa++) {
932 0 : if(pdb.getAtomName(frg_atoms[a])==his_n[aa]) {
933 0 : ri.atom[aa] = atm;
934 0 : break;
935 : }
936 : }
937 : }
938 0 : ri.numAtoms = 5;
939 0 : total_rings_atoms += 3;
940 0 : ri.rtype = RingInfo::R_HIS;
941 0 : ringInfo.push_back(ri);
942 : } else {
943 0 : plumed_merror("Unknown Ring Fragment: " + frg);
944 : }
945 : }
946 : }
947 :
948 10620 : for(unsigned cs=0; cs<chemicalshifts.size(); cs++) chemicalshifts[cs].csatoms += total_rings_atoms;
949 486 : }
950 :
951 18 : void CS2Backbone::calculate()
952 : {
953 18 : if(pbc) makeWhole();
954 18 : if(getExchangeStep()) box_count=0;
955 18 : if(box_count==0) update_neighb();
956 18 : compute_ring_parameters();
957 :
958 18 : std::vector<double> camshift_sigma2(6);
959 18 : camshift_sigma2[0] = 0.08; // HA
960 18 : camshift_sigma2[1] = 0.30; // HN
961 18 : camshift_sigma2[2] = 9.00; // NH
962 18 : camshift_sigma2[3] = 1.30; // CA
963 18 : camshift_sigma2[4] = 1.56; // CB
964 18 : camshift_sigma2[5] = 1.70; // CO
965 :
966 : std::vector<Vector> cs_derivs;
967 : std::vector<Vector> aa_derivs;
968 : std::vector<unsigned> cs_atoms;
969 : std::vector<double> all_shifts;
970 :
971 18 : cs_derivs.resize(chemicalshifts.size()*max_cs_atoms,Vector(0,0,0));
972 18 : cs_atoms.resize(chemicalshifts.size()*max_cs_atoms,0);
973 18 : all_shifts.resize(chemicalshifts.size(),0);
974 18 : if(camshift||getDoScore()) aa_derivs.resize(getNumberOfAtoms(),Vector(0,0,0));
975 :
976 18 : unsigned stride = comm.Get_size();
977 18 : unsigned rank = comm.Get_rank();
978 18 : if(serial) {
979 : stride = 1;
980 : rank = 0;
981 : }
982 :
983 18 : unsigned nt=OpenMP::getNumThreads();
984 18 : if(nt*stride*2>chemicalshifts.size()) nt=1;
985 :
986 : // a single loop over all chemical shifts
987 18 : #pragma omp parallel num_threads(nt)
988 : {
989 : #pragma omp for schedule(dynamic)
990 : for(unsigned cs=rank; cs<chemicalshifts.size(); cs+=stride) {
991 : const unsigned kdx=cs*max_cs_atoms;
992 : const ChemicalShift *myfrag = &chemicalshifts[cs];
993 : const unsigned aa_kind = myfrag->res_kind;
994 : const unsigned at_kind = myfrag->atm_kind;
995 :
996 : double shift = db.CONSTAAPREV(aa_kind,at_kind)[myfrag->res_type_prev] +
997 : db.CONSTAACURR(aa_kind,at_kind)[myfrag->res_type_curr] +
998 : db.CONSTAANEXT(aa_kind,at_kind)[myfrag->res_type_next];
999 :
1000 : const unsigned ipos = myfrag->ipos;
1001 : cs_atoms[kdx+0] = ipos;
1002 : unsigned atom_counter = 1;
1003 :
1004 : //BACKBONE (PREV CURR NEXT)
1005 : const double * CONST_BB2 = db.CONST_BB2(aa_kind,at_kind);
1006 : const unsigned bbsize = 16;
1007 : for(unsigned q=0; q<bbsize; q++) {
1008 : const double cb2q = CONST_BB2[q];
1009 : if(cb2q==0.) continue;
1010 : const unsigned jpos = myfrag->bb[q];
1011 : if(ipos==jpos) continue;
1012 : const Vector distance = delta(getPosition(jpos),getPosition(ipos));
1013 : const double d = distance.modulo();
1014 : const double fact = cb2q/d;
1015 :
1016 : shift += cb2q*d;
1017 : const Vector der = fact*distance;
1018 :
1019 : cs_derivs[kdx+0] += der;
1020 : cs_derivs[kdx+q+atom_counter] = -der;
1021 : cs_atoms[kdx+q+atom_counter] = jpos;
1022 : }
1023 :
1024 : atom_counter += bbsize;
1025 :
1026 : //DIHEDRAL ANGLES
1027 : const double *CO_DA = db.CO_DA(aa_kind,at_kind);
1028 : //Phi
1029 : {
1030 : const Vector d0 = delta(getPosition(myfrag->bb[Nc]), getPosition(myfrag->bb[Cp]));
1031 : const Vector d1 = delta(getPosition(myfrag->bb[CAc]), getPosition(myfrag->bb[Nc]));
1032 : const Vector d2 = delta(getPosition(myfrag->bb[Cc]), getPosition(myfrag->bb[CAc]));
1033 : Torsion t;
1034 : Vector dd0, dd1, dd2;
1035 : const double t_phi = t.compute(d0,d1,d2,dd0,dd1,dd2);
1036 : const double *PARS_DA = db.PARS_DA(aa_kind,at_kind,0);
1037 : const double val1 = 3.*t_phi+PARS_DA[3];
1038 : const double val2 = t_phi+PARS_DA[4];
1039 : shift += CO_DA[0]*(PARS_DA[0]*std::cos(val1)+PARS_DA[1]*std::cos(val2)+PARS_DA[2]);
1040 : const double fact = -CO_DA[0]*(+3.*PARS_DA[0]*std::sin(val1)+PARS_DA[1]*std::sin(val2));
1041 :
1042 : cs_derivs[kdx+Cp+1] += fact*dd0;
1043 : cs_derivs[kdx+Nc+1] += fact*(dd1-dd0);
1044 : cs_derivs[kdx+CAc+1]+= fact*(dd2-dd1);
1045 : cs_derivs[kdx+Cc+1] += -fact*dd2;
1046 : cs_atoms[kdx+Cp+1] = myfrag->bb[Cp];
1047 : cs_atoms[kdx+Nc+1] = myfrag->bb[Nc];
1048 : cs_atoms[kdx+CAc+1]= myfrag->bb[CAc];
1049 : cs_atoms[kdx+Cc+1] = myfrag->bb[Cc];
1050 : }
1051 :
1052 : //Psi
1053 : {
1054 : const Vector d0 = delta(getPosition(myfrag->bb[CAc]), getPosition(myfrag->bb[Nc]));
1055 : const Vector d1 = delta(getPosition(myfrag->bb[Cc]), getPosition(myfrag->bb[CAc]));
1056 : const Vector d2 = delta(getPosition(myfrag->bb[Nn]), getPosition(myfrag->bb[Cc]));
1057 : Torsion t;
1058 : Vector dd0, dd1, dd2;
1059 : const double t_psi = t.compute(d0,d1,d2,dd0,dd1,dd2);
1060 : const double *PARS_DA = db.PARS_DA(aa_kind,at_kind,1);
1061 : const double val1 = 3.*t_psi+PARS_DA[3];
1062 : const double val2 = t_psi+PARS_DA[4];
1063 : shift += CO_DA[1]*(PARS_DA[0]*std::cos(val1)+PARS_DA[1]*std::cos(val2)+PARS_DA[2]);
1064 : const double fact = -CO_DA[1]*(+3.*PARS_DA[0]*std::sin(val1)+PARS_DA[1]*std::sin(val2));
1065 :
1066 : cs_derivs[kdx+Nc+1] += fact*dd0;
1067 : cs_derivs[kdx+CAc+1] += fact*(dd1-dd0);
1068 : cs_derivs[kdx+Cc+1] += fact*(dd2-dd1);
1069 : cs_derivs[kdx+Nn+1] += -fact*dd2;
1070 : cs_atoms[kdx+Nc+1] = myfrag->bb[Nc];
1071 : cs_atoms[kdx+CAc+1]= myfrag->bb[CAc];
1072 : cs_atoms[kdx+Cc+1] = myfrag->bb[Cc];
1073 : cs_atoms[kdx+Nn+1] = myfrag->bb[Nn];
1074 : }
1075 :
1076 : //Chi
1077 : if(myfrag->has_chi1) {
1078 : const Vector d0 = delta(getPosition(myfrag->bb[CAc]), getPosition(myfrag->bb[Nc]));
1079 : const Vector d1 = delta(getPosition(myfrag->bb[CBc]), getPosition(myfrag->bb[CAc]));
1080 : const Vector d2 = delta(getPosition(myfrag->bb[CGc]), getPosition(myfrag->bb[CBc]));
1081 : Torsion t;
1082 : Vector dd0, dd1, dd2;
1083 : const double t_chi1 = t.compute(d0,d1,d2,dd0,dd1,dd2);
1084 : const double *PARS_DA = db.PARS_DA(aa_kind,at_kind,2);
1085 : const double val1 = 3.*t_chi1+PARS_DA[3];
1086 : const double val2 = t_chi1+PARS_DA[4];
1087 : shift += CO_DA[2]*(PARS_DA[0]*std::cos(val1)+PARS_DA[1]*std::cos(val2)+PARS_DA[2]);
1088 : const double fact = -CO_DA[2]*(+3.*PARS_DA[0]*std::sin(val1)+PARS_DA[1]*std::sin(val2));
1089 :
1090 : cs_derivs[kdx+Nc+1] += fact*dd0;
1091 : cs_derivs[kdx+CAc+1] += fact*(dd1-dd0);
1092 : cs_derivs[kdx+CBc+1] += fact*(dd2-dd1);
1093 : cs_derivs[kdx+CGc+1] += -fact*dd2;
1094 : cs_atoms[kdx+Nc+1] = myfrag->bb[Nc];
1095 : cs_atoms[kdx+CAc+1] = myfrag->bb[CAc];
1096 : cs_atoms[kdx+CBc+1] = myfrag->bb[CBc];
1097 : cs_atoms[kdx+CGc+1] = myfrag->bb[CGc];
1098 :
1099 : atom_counter += 2;
1100 : }
1101 : //END OF DIHE
1102 :
1103 : //SIDE CHAIN
1104 : const double * CONST_SC2 = db.CONST_SC2(aa_kind,at_kind,myfrag->res_type_curr);
1105 : const unsigned sidsize = myfrag->side_chain.size();
1106 : for(unsigned q=0; q<sidsize; q++) {
1107 : const double cs2q = CONST_SC2[q];
1108 : if(cs2q==0.) continue;
1109 : const unsigned jpos = myfrag->side_chain[q];
1110 : if(ipos==jpos) continue;
1111 : const Vector distance = delta(getPosition(jpos),getPosition(ipos));
1112 : const double d = distance.modulo();
1113 : const double fact = cs2q/d;
1114 :
1115 : shift += cs2q*d;
1116 : const Vector der = fact*distance;
1117 : cs_derivs[kdx+0] += der;
1118 : cs_derivs[kdx+q+atom_counter] = -der;
1119 : cs_atoms[kdx+q+atom_counter] = jpos;
1120 : }
1121 :
1122 : atom_counter += sidsize;
1123 :
1124 : //EXTRA DIST
1125 : const double * CONST_XD = db.CONST_XD(aa_kind,at_kind);
1126 : const unsigned xdsize=myfrag->xd1.size();
1127 : for(unsigned q=0; q<xdsize; q++) {
1128 : const double cxdq = CONST_XD[q];
1129 : if(cxdq==0.) continue;
1130 : if(myfrag->xd1[q]==-1||myfrag->xd2[q]==-1) continue;
1131 : const Vector distance = delta(getPosition(myfrag->xd1[q]),getPosition(myfrag->xd2[q]));
1132 : const double d = distance.modulo();
1133 : const double fact = cxdq/d;
1134 :
1135 : shift += cxdq*d;
1136 : const Vector der = fact*distance;
1137 : cs_derivs[kdx+2*q+atom_counter ] = der;
1138 : cs_derivs[kdx+2*q+atom_counter+1] = -der;
1139 : cs_atoms[kdx+2*q+atom_counter] = myfrag->xd2[q];
1140 : cs_atoms[kdx+2*q+atom_counter+1] = myfrag->xd1[q];
1141 : }
1142 :
1143 : atom_counter += 2*xdsize;
1144 :
1145 : //RINGS
1146 : const double *rc = db.CO_RING(aa_kind,at_kind);
1147 : const unsigned rsize = ringInfo.size();
1148 : // cycle over the list of rings
1149 : for(unsigned q=0; q<rsize; q++) {
1150 : // compute angle from ring middle point to current atom position
1151 : // get distance std::vector from query atom to ring center and normal std::vector to ring plane
1152 : const Vector n = ringInfo[q].normVect;
1153 : const double nL = ringInfo[q].lengthNV;
1154 : const double inL2 = ringInfo[q].lengthN2;
1155 :
1156 : const Vector d = delta(ringInfo[q].position, getPosition(ipos));
1157 : const double dL2 = d.modulo2();
1158 : double dL = std::sqrt(dL2);
1159 : const double idL3 = 1./(dL2*dL);
1160 :
1161 : const double dn = dotProduct(d,n);
1162 : const double dn2 = dn*dn;
1163 : const double dLnL = dL*nL;
1164 : const double dL_nL = dL/nL;
1165 :
1166 : const double ang2 = dn2*inL2/dL2;
1167 : const double u = 1.-3.*ang2;
1168 : const double cc = rc[ringInfo[q].rtype];
1169 :
1170 : shift += cc*u*idL3;
1171 :
1172 : const double fUU = -6.*dn*inL2;
1173 : const double fUQ = fUU/dL;
1174 : const Vector gradUQ = fUQ*(dL2*n - dn*d);
1175 : const Vector gradVQ = (3.*dL*u)*d;
1176 :
1177 : const double fact = cc*idL3*idL3;
1178 : cs_derivs[kdx+0] += fact*(gradUQ - gradVQ);
1179 :
1180 : const double fU = fUU/nL;
1181 : double OneOverN = 1./6.;
1182 : if(ringInfo[q].numAtoms==5) OneOverN=1./3.;
1183 : const Vector factor2 = OneOverN*n;
1184 : const Vector factor4 = (OneOverN/dL_nL)*d;
1185 :
1186 : const Vector gradV = -OneOverN*gradVQ;
1187 :
1188 : if(ringInfo[q].numAtoms==6) {
1189 : // update forces on ring atoms
1190 : for(unsigned at=0; at<6; at++) {
1191 : const Vector ab = crossProduct(d,ringInfo[q].g[at]);
1192 : const Vector c = crossProduct(n,ringInfo[q].g[at]);
1193 : const Vector factor3 = 0.5*dL_nL*c;
1194 : const Vector factor1 = 0.5*ab;
1195 : const Vector gradU = fU*( dLnL*(factor1 - factor2) -dn*(factor3 - factor4) );
1196 : cs_derivs[kdx+at+atom_counter] = fact*(gradU - gradV);
1197 : cs_atoms[kdx+at+atom_counter] = ringInfo[q].atom[at];
1198 : }
1199 : atom_counter += 6;
1200 : } else {
1201 : for(unsigned at=0; at<3; at++) {
1202 : const Vector ab = crossProduct(d,ringInfo[q].g[at]);
1203 : const Vector c = crossProduct(n,ringInfo[q].g[at]);
1204 : const Vector factor3 = dL_nL*c;
1205 : const Vector factor1 = ab;
1206 : const Vector gradU = fU*( dLnL*(factor1 - factor2) -dn*(factor3 - factor4) );
1207 : cs_derivs[kdx+at+atom_counter] = fact*(gradU - gradV);
1208 : }
1209 : cs_atoms[kdx+atom_counter] = ringInfo[q].atom[0];
1210 : cs_atoms[kdx+atom_counter+1] = ringInfo[q].atom[2];
1211 : cs_atoms[kdx+atom_counter+2] = ringInfo[q].atom[3];
1212 : atom_counter += 3;
1213 : }
1214 : }
1215 : //END OF RINGS
1216 :
1217 : //NON BOND
1218 : const double * CONST_CO_SPHERE3 = db.CO_SPHERE(aa_kind,at_kind,0);
1219 : const double * CONST_CO_SPHERE = db.CO_SPHERE(aa_kind,at_kind,1);
1220 : const unsigned boxsize = myfrag->box_nb.size();
1221 : for(unsigned q=0; q<boxsize; q++) {
1222 : const unsigned jpos = myfrag->box_nb[q];
1223 : const Vector distance = delta(getPosition(jpos),getPosition(ipos));
1224 : const double d2 = distance.modulo2();
1225 :
1226 : if(d2<cutOffDist2) {
1227 : double factor1 = std::sqrt(d2);
1228 : double dfactor1 = 1./factor1;
1229 : double factor3 = dfactor1*dfactor1*dfactor1;
1230 : double dfactor3 = -3.*factor3*dfactor1*dfactor1;
1231 :
1232 : if(d2>cutOnDist2) {
1233 : const double af = cutOffDist2 - d2;
1234 : const double bf = cutOffDist2 - 3.*cutOnDist2 + 2.*d2;
1235 : const double cf = invswitch*af;
1236 : const double df = cf*af*bf;
1237 : factor1 *= df;
1238 : factor3 *= df;
1239 :
1240 : const double d4 = d2*d2;
1241 : const double af1 = 15.*cutOnDist2*d2;
1242 : const double bf1 = -14.*d4;
1243 : const double cf1 = -3.*cutOffDist2*cutOnDist2 + cutOffDist2*d2;
1244 : const double df1 = af1+bf1+cf1;
1245 : dfactor1 *= cf*(cutOffDist4+df1);
1246 :
1247 : const double af3 = +2.*cutOffDist2*cutOnDist2;
1248 : const double bf3 = d2*(cutOffDist2+cutOnDist2);
1249 : const double cf3 = -2.*d4;
1250 : const double df3 = (af3+bf3+cf3)*d2;
1251 : dfactor3 *= invswitch*(cutMixed+df3);
1252 : }
1253 :
1254 : const unsigned t = type[jpos];
1255 : shift += factor1*CONST_CO_SPHERE[t] + factor3*CONST_CO_SPHERE3[t] ;
1256 : const double fact = dfactor1*CONST_CO_SPHERE[t]+dfactor3*CONST_CO_SPHERE3[t];
1257 : const Vector der = fact*distance;
1258 :
1259 : cs_derivs[kdx+0] += der;
1260 : cs_derivs[kdx+q+atom_counter] = -der;
1261 : cs_atoms[kdx+q+atom_counter] = jpos;
1262 : }
1263 : }
1264 : //END NON BOND
1265 :
1266 : atom_counter += boxsize;
1267 : all_shifts[cs] = shift;
1268 : }
1269 : }
1270 :
1271 18 : ++box_count;
1272 18 : if(box_count == box_nupdate) box_count = 0;
1273 :
1274 18 : if(!camshift) {
1275 13 : if(!serial) {
1276 13 : if(!getDoScore()) {
1277 9 : comm.Sum(&cs_derivs[0][0], 3*cs_derivs.size());
1278 9 : comm.Sum(&cs_atoms[0], cs_atoms.size());
1279 : }
1280 13 : comm.Sum(&all_shifts[0], chemicalshifts.size());
1281 : }
1282 7670 : for(unsigned cs=0; cs<chemicalshifts.size(); cs++) {
1283 7657 : Value *comp = chemicalshifts[cs].comp;
1284 7657 : comp->set(all_shifts[cs]);
1285 7657 : if(getDoScore()) setCalcData(cs, all_shifts[cs]);
1286 : else {
1287 5301 : const unsigned kdx=cs*max_cs_atoms;
1288 5301 : Tensor csvirial;
1289 1270127 : for(unsigned i=0; i<chemicalshifts[cs].totcsatoms; i++) {
1290 1264826 : setAtomsDerivatives(comp,cs_atoms[kdx+i],cs_derivs[kdx+i]);
1291 1264826 : csvirial-=Tensor(getPosition(cs_atoms[kdx+i]),cs_derivs[kdx+i]);
1292 : }
1293 5301 : setBoxDerivatives(comp,csvirial);
1294 : }
1295 : }
1296 13 : if(!getDoScore()) return;
1297 : }
1298 :
1299 9 : double score = 0.;
1300 :
1301 : /* Metainference */
1302 9 : if(getDoScore()) {
1303 4 : score = getScore();
1304 1182 : for(unsigned cs=rank; cs<chemicalshifts.size(); cs+=stride) {
1305 1178 : const unsigned kdx=cs*max_cs_atoms;
1306 : const double fact = getMetaDer(cs);
1307 282215 : for(unsigned i=0; i<chemicalshifts[cs].totcsatoms; i++) {
1308 281037 : aa_derivs[cs_atoms[kdx+i]] += cs_derivs[kdx+i]*fact;
1309 : }
1310 : }
1311 : }
1312 :
1313 : /* camshift */
1314 9 : if(camshift) {
1315 1772 : for(unsigned cs=rank; cs<chemicalshifts.size(); cs+=stride) {
1316 1767 : const unsigned kdx=cs*max_cs_atoms;
1317 1767 : score += (all_shifts[cs] - chemicalshifts[cs].exp_cs)*(all_shifts[cs] - chemicalshifts[cs].exp_cs)/camshift_sigma2[chemicalshifts[cs].atm_kind];
1318 1767 : double fact = 2.0*(all_shifts[cs] - chemicalshifts[cs].exp_cs)/camshift_sigma2[chemicalshifts[cs].atm_kind];
1319 423482 : for(unsigned i=0; i<chemicalshifts[cs].totcsatoms; i++) {
1320 421715 : aa_derivs[cs_atoms[kdx+i]] += cs_derivs[kdx+i]*fact;
1321 : }
1322 : }
1323 : }
1324 :
1325 9 : if(!serial) {
1326 9 : comm.Sum(&aa_derivs[0][0], 3*aa_derivs.size());
1327 9 : if(camshift) comm.Sum(&score, 1);
1328 : }
1329 :
1330 9 : Tensor virial;
1331 13069 : for(unsigned i=rank; i<getNumberOfAtoms(); i+=stride) {
1332 13060 : virial += Tensor(getPosition(i), aa_derivs[i]);
1333 : }
1334 :
1335 9 : if(!serial) {
1336 9 : comm.Sum(&virial[0][0], 9);
1337 : }
1338 :
1339 : /* calculate final derivatives */
1340 : Value* val;
1341 9 : if(getDoScore()) {
1342 4 : val=getPntrToComponent("score");
1343 4 : setScore(score);
1344 : } else {
1345 5 : val=getPntrToValue();
1346 5 : setValue(score);
1347 : }
1348 :
1349 : /* at this point we cycle over all atoms */
1350 23517 : for(unsigned i=0; i<getNumberOfAtoms(); i++) setAtomsDerivatives(val, i, aa_derivs[i]);
1351 9 : setBoxDerivatives(val,-virial);
1352 : }
1353 :
1354 18 : void CS2Backbone::update_neighb() {
1355 : // cycle over chemical shifts
1356 18 : unsigned nt=OpenMP::getNumThreads();
1357 18 : #pragma omp parallel for num_threads(nt)
1358 : for(unsigned cs=0; cs<chemicalshifts.size(); cs++) {
1359 : const unsigned boxsize = getNumberOfAtoms();
1360 : chemicalshifts[cs].box_nb.clear();
1361 : chemicalshifts[cs].box_nb.reserve(150);
1362 : const unsigned res_curr = res_num[chemicalshifts[cs].ipos];
1363 : for(unsigned bat=0; bat<boxsize; bat++) {
1364 : const unsigned res_dist = std::abs(static_cast<int>(res_curr-res_num[bat]));
1365 : if(res_dist<2) continue;
1366 : const Vector distance = delta(getPosition(bat),getPosition(chemicalshifts[cs].ipos));
1367 : const double d2=distance.modulo2();
1368 : if(d2<cutOffNB2) chemicalshifts[cs].box_nb.push_back(bat);
1369 : }
1370 : chemicalshifts[cs].totcsatoms = chemicalshifts[cs].csatoms + chemicalshifts[cs].box_nb.size();
1371 : }
1372 18 : max_cs_atoms=0;
1373 10620 : for(unsigned cs=0; cs<chemicalshifts.size(); cs++) {
1374 10602 : if(chemicalshifts[cs].totcsatoms>max_cs_atoms) max_cs_atoms = chemicalshifts[cs].totcsatoms;
1375 : }
1376 18 : }
1377 :
1378 18 : void CS2Backbone::compute_ring_parameters() {
1379 378 : for(unsigned i=0; i<ringInfo.size(); i++) {
1380 360 : const unsigned size = ringInfo[i].numAtoms;
1381 360 : if(size==6) {
1382 342 : ringInfo[i].g[0] = delta(getPosition(ringInfo[i].atom[4]),getPosition(ringInfo[i].atom[2]));
1383 342 : ringInfo[i].g[1] = delta(getPosition(ringInfo[i].atom[5]),getPosition(ringInfo[i].atom[3]));
1384 342 : ringInfo[i].g[2] = delta(getPosition(ringInfo[i].atom[0]),getPosition(ringInfo[i].atom[4]));
1385 342 : ringInfo[i].g[3] = delta(getPosition(ringInfo[i].atom[1]),getPosition(ringInfo[i].atom[5]));
1386 342 : ringInfo[i].g[4] = delta(getPosition(ringInfo[i].atom[2]),getPosition(ringInfo[i].atom[0]));
1387 342 : ringInfo[i].g[5] = delta(getPosition(ringInfo[i].atom[3]),getPosition(ringInfo[i].atom[1]));
1388 : // ring center
1389 342 : Vector midP = getPosition(ringInfo[i].atom[0]);
1390 2052 : for(unsigned j=1; j<size; j++) midP += getPosition(ringInfo[i].atom[j]);
1391 342 : ringInfo[i].position = midP/6.;
1392 : // compute normal std::vector to plane
1393 342 : Vector n1 = crossProduct(ringInfo[i].g[2], -ringInfo[i].g[4]);
1394 342 : Vector n2 = crossProduct(ringInfo[i].g[5], -ringInfo[i].g[1]);
1395 342 : ringInfo[i].normVect = 0.5*(n1 + n2);
1396 : } else {
1397 18 : ringInfo[i].g[0] = delta(getPosition(ringInfo[i].atom[3]),getPosition(ringInfo[i].atom[2]));
1398 18 : ringInfo[i].g[1] = delta(getPosition(ringInfo[i].atom[0]),getPosition(ringInfo[i].atom[3]));
1399 18 : ringInfo[i].g[2] = delta(getPosition(ringInfo[i].atom[2]),getPosition(ringInfo[i].atom[0]));
1400 : // ring center
1401 18 : ringInfo[i].position = (getPosition(ringInfo[i].atom[0])+getPosition(ringInfo[i].atom[2])+getPosition(ringInfo[i].atom[3]))/3.;
1402 : // ring plane normal std::vector
1403 18 : ringInfo[i].normVect = crossProduct(ringInfo[i].g[1],-ringInfo[i].g[2]);
1404 :
1405 : }
1406 : // calculate squared length and length of normal std::vector
1407 360 : ringInfo[i].lengthN2 = 1./ringInfo[i].normVect.modulo2();
1408 360 : ringInfo[i].lengthNV = 1./std::sqrt(ringInfo[i].lengthN2);
1409 : }
1410 18 : }
1411 :
1412 31806 : CS2Backbone::aa_t CS2Backbone::frag2enum(const std::string &aa) {
1413 : aa_t type = ALA;
1414 31806 : if (aa == "ALA") type = ALA;
1415 29934 : else if (aa == "ARG") type = ARG;
1416 28548 : else if (aa == "ASN") type = ASN;
1417 26910 : else if (aa == "ASP") type = ASP;
1418 25380 : else if (aa == "ASH") type = ASP;
1419 25380 : else if (aa == "CYS") type = CYS;
1420 24858 : else if (aa == "CYM") type = CYS;
1421 24858 : else if (aa == "GLN") type = GLN;
1422 24390 : else if (aa == "GLU") type = GLU;
1423 22068 : else if (aa == "GLH") type = GLU;
1424 22068 : else if (aa == "GLY") type = GLY;
1425 16974 : else if (aa == "HIS") type = HIS;
1426 16974 : else if (aa == "HSE") type = HIS;
1427 16974 : else if (aa == "HIE") type = HIS;
1428 16974 : else if (aa == "HSP") type = HIS;
1429 16974 : else if (aa == "HIP") type = HIS;
1430 16974 : else if (aa == "HSD") type = HIS;
1431 16974 : else if (aa == "HID") type = HIS;
1432 16974 : else if (aa == "ILE") type = ILE;
1433 15174 : else if (aa == "LEU") type = LEU;
1434 13536 : else if (aa == "LYS") type = LYS;
1435 10746 : else if (aa == "MET") type = MET;
1436 9936 : else if (aa == "PHE") type = PHE;
1437 6876 : else if (aa == "PRO") type = PRO;
1438 5940 : else if (aa == "SER") type = SER;
1439 4302 : else if (aa == "THR") type = THR;
1440 2196 : else if (aa == "TRP") type = TRP;
1441 1980 : else if (aa == "TYR") type = TYR;
1442 1602 : else if (aa == "VAL") type = VAL;
1443 0 : else if (aa == "UNK") type = UNK;
1444 0 : else plumed_merror("Error converting std::string " + aa + " into amino acid index: not a valid 3-letter code");
1445 31806 : return type;
1446 : }
1447 :
1448 10602 : std::vector<std::string> CS2Backbone::side_chain_atoms(const std::string &s) {
1449 : std::vector<std::string> sc;
1450 :
1451 10602 : if(s=="ALA") {
1452 648 : sc.push_back( "CB" );
1453 648 : sc.push_back( "HB1" );
1454 648 : sc.push_back( "HB2" );
1455 648 : sc.push_back( "HB3" );
1456 648 : return sc;
1457 9954 : } else if(s=="ARG") {
1458 468 : sc.push_back( "CB" );
1459 468 : sc.push_back( "CG" );
1460 468 : sc.push_back( "CD" );
1461 468 : sc.push_back( "NE" );
1462 468 : sc.push_back( "CZ" );
1463 468 : sc.push_back( "NH1" );
1464 468 : sc.push_back( "NH2" );
1465 468 : sc.push_back( "NH3" );
1466 468 : sc.push_back( "HB1" );
1467 468 : sc.push_back( "HB2" );
1468 468 : sc.push_back( "HB3" );
1469 468 : sc.push_back( "HG1" );
1470 468 : sc.push_back( "HG2" );
1471 468 : sc.push_back( "HG3" );
1472 468 : sc.push_back( "HD1" );
1473 468 : sc.push_back( "HD2" );
1474 468 : sc.push_back( "HD3" );
1475 468 : sc.push_back( "HE" );
1476 468 : sc.push_back( "HH11" );
1477 468 : sc.push_back( "HH12" );
1478 468 : sc.push_back( "HH21" );
1479 468 : sc.push_back( "HH22" );
1480 468 : sc.push_back( "1HH1" );
1481 468 : sc.push_back( "2HH1" );
1482 468 : sc.push_back( "1HH2" );
1483 468 : sc.push_back( "2HH2" );
1484 468 : return sc;
1485 9486 : } else if(s=="ASN") {
1486 594 : sc.push_back( "CB" );
1487 594 : sc.push_back( "CG" );
1488 594 : sc.push_back( "OD1" );
1489 594 : sc.push_back( "ND2" );
1490 594 : sc.push_back( "HB1" );
1491 594 : sc.push_back( "HB2" );
1492 594 : sc.push_back( "HB3" );
1493 594 : sc.push_back( "HD21" );
1494 594 : sc.push_back( "HD22" );
1495 594 : sc.push_back( "1HD2" );
1496 594 : sc.push_back( "2HD2" );
1497 594 : return sc;
1498 17226 : } else if(s=="ASP"||s=="ASH") {
1499 558 : sc.push_back( "CB" );
1500 558 : sc.push_back( "CG" );
1501 558 : sc.push_back( "OD1" );
1502 558 : sc.push_back( "OD2" );
1503 558 : sc.push_back( "HB1" );
1504 558 : sc.push_back( "HB2" );
1505 558 : sc.push_back( "HB3" );
1506 558 : return sc;
1507 16668 : } else if(s=="CYS"||s=="CYM") {
1508 0 : sc.push_back( "CB" );
1509 0 : sc.push_back( "SG" );
1510 0 : sc.push_back( "HB1" );
1511 0 : sc.push_back( "HB2" );
1512 0 : sc.push_back( "HB3" );
1513 0 : sc.push_back( "HG1" );
1514 0 : sc.push_back( "HG" );
1515 0 : return sc;
1516 8334 : } else if(s=="GLN") {
1517 162 : sc.push_back( "CB" );
1518 162 : sc.push_back( "CG" );
1519 162 : sc.push_back( "CD" );
1520 162 : sc.push_back( "OE1" );
1521 162 : sc.push_back( "NE2" );
1522 162 : sc.push_back( "HB1" );
1523 162 : sc.push_back( "HB2" );
1524 162 : sc.push_back( "HB3" );
1525 162 : sc.push_back( "HG1" );
1526 162 : sc.push_back( "HG2" );
1527 162 : sc.push_back( "HG3" );
1528 162 : sc.push_back( "HE21" );
1529 162 : sc.push_back( "HE22" );
1530 162 : sc.push_back( "1HE2" );
1531 162 : sc.push_back( "2HE2" );
1532 162 : return sc;
1533 15552 : } else if(s=="GLU"||s=="GLH") {
1534 792 : sc.push_back( "CB" );
1535 792 : sc.push_back( "CG" );
1536 792 : sc.push_back( "CD" );
1537 792 : sc.push_back( "OE1" );
1538 792 : sc.push_back( "OE2" );
1539 792 : sc.push_back( "HB1" );
1540 792 : sc.push_back( "HB2" );
1541 792 : sc.push_back( "HB3" );
1542 792 : sc.push_back( "HG1" );
1543 792 : sc.push_back( "HG2" );
1544 792 : sc.push_back( "HG3" );
1545 792 : return sc;
1546 7380 : } else if(s=="GLY") {
1547 1494 : sc.push_back( "HA2" );
1548 1494 : return sc;
1549 41202 : } else if(s=="HIS"||s=="HSE"||s=="HIE"||s=="HSD"||s=="HID"||s=="HIP"||s=="HSP") {
1550 0 : sc.push_back( "CB" );
1551 0 : sc.push_back( "CG" );
1552 0 : sc.push_back( "ND1" );
1553 0 : sc.push_back( "CD2" );
1554 0 : sc.push_back( "CE1" );
1555 0 : sc.push_back( "NE2" );
1556 0 : sc.push_back( "HB1" );
1557 0 : sc.push_back( "HB2" );
1558 0 : sc.push_back( "HB3" );
1559 0 : sc.push_back( "HD1" );
1560 0 : sc.push_back( "HD2" );
1561 0 : sc.push_back( "HE1" );
1562 0 : sc.push_back( "HE2" );
1563 0 : return sc;
1564 5886 : } else if(s=="ILE") {
1565 540 : sc.push_back( "CB" );
1566 540 : sc.push_back( "CG1" );
1567 540 : sc.push_back( "CG2" );
1568 540 : sc.push_back( "CD" );
1569 540 : sc.push_back( "HB" );
1570 540 : sc.push_back( "HG11" );
1571 540 : sc.push_back( "HG12" );
1572 540 : sc.push_back( "HG21" );
1573 540 : sc.push_back( "HG22" );
1574 540 : sc.push_back( "HG23" );
1575 540 : sc.push_back( "1HG1" );
1576 540 : sc.push_back( "2HG1" );
1577 540 : sc.push_back( "1HG2" );
1578 540 : sc.push_back( "2HG2" );
1579 540 : sc.push_back( "3HG2" );
1580 540 : sc.push_back( "HD1" );
1581 540 : sc.push_back( "HD2" );
1582 540 : sc.push_back( "HD3" );
1583 540 : return sc;
1584 5346 : } else if(s=="LEU") {
1585 648 : sc.push_back( "CB" );
1586 648 : sc.push_back( "CG" );
1587 648 : sc.push_back( "CD1" );
1588 648 : sc.push_back( "CD2" );
1589 648 : sc.push_back( "HB1" );
1590 648 : sc.push_back( "HB2" );
1591 648 : sc.push_back( "HB3" );
1592 648 : sc.push_back( "HG" );
1593 648 : sc.push_back( "HD11" );
1594 648 : sc.push_back( "HD12" );
1595 648 : sc.push_back( "HD13" );
1596 648 : sc.push_back( "HD21" );
1597 648 : sc.push_back( "HD22" );
1598 648 : sc.push_back( "HD23" );
1599 648 : sc.push_back( "1HD1" );
1600 648 : sc.push_back( "2HD1" );
1601 648 : sc.push_back( "3HD1" );
1602 648 : sc.push_back( "1HD2" );
1603 648 : sc.push_back( "2HD2" );
1604 648 : sc.push_back( "3HD2" );
1605 648 : return sc;
1606 4698 : } else if(s=="LYS") {
1607 1008 : sc.push_back( "CB" );
1608 1008 : sc.push_back( "CG" );
1609 1008 : sc.push_back( "CD" );
1610 1008 : sc.push_back( "CE" );
1611 1008 : sc.push_back( "NZ" );
1612 1008 : sc.push_back( "HB1" );
1613 1008 : sc.push_back( "HB2" );
1614 1008 : sc.push_back( "HB3" );
1615 1008 : sc.push_back( "HG1" );
1616 1008 : sc.push_back( "HG2" );
1617 1008 : sc.push_back( "HG3" );
1618 1008 : sc.push_back( "HD1" );
1619 1008 : sc.push_back( "HD2" );
1620 1008 : sc.push_back( "HD3" );
1621 1008 : sc.push_back( "HE1" );
1622 1008 : sc.push_back( "HE2" );
1623 1008 : sc.push_back( "HE3" );
1624 1008 : sc.push_back( "HZ1" );
1625 1008 : sc.push_back( "HZ2" );
1626 1008 : sc.push_back( "HZ3" );
1627 1008 : return sc;
1628 3690 : } else if(s=="MET") {
1629 288 : sc.push_back( "CB" );
1630 288 : sc.push_back( "CG" );
1631 288 : sc.push_back( "SD" );
1632 288 : sc.push_back( "CE" );
1633 288 : sc.push_back( "HB1" );
1634 288 : sc.push_back( "HB2" );
1635 288 : sc.push_back( "HB3" );
1636 288 : sc.push_back( "HG1" );
1637 288 : sc.push_back( "HG2" );
1638 288 : sc.push_back( "HG3" );
1639 288 : sc.push_back( "HE1" );
1640 288 : sc.push_back( "HE2" );
1641 288 : sc.push_back( "HE3" );
1642 288 : return sc;
1643 3402 : } else if(s=="PHE") {
1644 1098 : sc.push_back( "CB" );
1645 1098 : sc.push_back( "CG" );
1646 1098 : sc.push_back( "CD1" );
1647 1098 : sc.push_back( "CD2" );
1648 1098 : sc.push_back( "CE1" );
1649 1098 : sc.push_back( "CE2" );
1650 1098 : sc.push_back( "CZ" );
1651 1098 : sc.push_back( "HB1" );
1652 1098 : sc.push_back( "HB2" );
1653 1098 : sc.push_back( "HB3" );
1654 1098 : sc.push_back( "HD1" );
1655 1098 : sc.push_back( "HD2" );
1656 1098 : sc.push_back( "HD3" );
1657 1098 : sc.push_back( "HE1" );
1658 1098 : sc.push_back( "HE2" );
1659 1098 : sc.push_back( "HE3" );
1660 1098 : sc.push_back( "HZ" );
1661 1098 : return sc;
1662 2304 : } else if(s=="PRO") {
1663 108 : sc.push_back( "CB" );
1664 108 : sc.push_back( "CG" );
1665 108 : sc.push_back( "CD" );
1666 108 : sc.push_back( "HB1" );
1667 108 : sc.push_back( "HB2" );
1668 108 : sc.push_back( "HB3" );
1669 108 : sc.push_back( "HG1" );
1670 108 : sc.push_back( "HG2" );
1671 108 : sc.push_back( "HG3" );
1672 108 : sc.push_back( "HD1" );
1673 108 : sc.push_back( "HD2" );
1674 108 : sc.push_back( "HD3" );
1675 108 : return sc;
1676 2196 : } else if(s=="SER") {
1677 630 : sc.push_back( "CB" );
1678 630 : sc.push_back( "OG" );
1679 630 : sc.push_back( "HB1" );
1680 630 : sc.push_back( "HB2" );
1681 630 : sc.push_back( "HB3" );
1682 630 : sc.push_back( "HG1" );
1683 630 : sc.push_back( "HG" );
1684 630 : return sc;
1685 1566 : } else if(s=="THR") {
1686 774 : sc.push_back( "CB" );
1687 774 : sc.push_back( "OG1" );
1688 774 : sc.push_back( "CG2" );
1689 774 : sc.push_back( "HB" );
1690 774 : sc.push_back( "HG1" );
1691 774 : sc.push_back( "HG21" );
1692 774 : sc.push_back( "HG22" );
1693 774 : sc.push_back( "HG23" );
1694 774 : sc.push_back( "1HG2" );
1695 774 : sc.push_back( "2HG2" );
1696 774 : sc.push_back( "3HG2" );
1697 774 : return sc;
1698 792 : } else if(s=="TRP") {
1699 72 : sc.push_back( "CB" );
1700 72 : sc.push_back( "CG" );
1701 72 : sc.push_back( "CD1" );
1702 72 : sc.push_back( "CD2" );
1703 72 : sc.push_back( "NE1" );
1704 72 : sc.push_back( "CE2" );
1705 72 : sc.push_back( "CE3" );
1706 72 : sc.push_back( "CZ2" );
1707 72 : sc.push_back( "CZ3" );
1708 72 : sc.push_back( "CH2" );
1709 72 : sc.push_back( "HB1" );
1710 72 : sc.push_back( "HB2" );
1711 72 : sc.push_back( "HB3" );
1712 72 : sc.push_back( "HD1" );
1713 72 : sc.push_back( "HE1" );
1714 72 : sc.push_back( "HE3" );
1715 72 : sc.push_back( "HZ2" );
1716 72 : sc.push_back( "HZ3" );
1717 72 : sc.push_back( "HH2" );
1718 72 : return sc;
1719 720 : } else if(s=="TYR") {
1720 144 : sc.push_back( "CB" );
1721 144 : sc.push_back( "CG" );
1722 144 : sc.push_back( "CD1" );
1723 144 : sc.push_back( "CD2" );
1724 144 : sc.push_back( "CE1" );
1725 144 : sc.push_back( "CE2" );
1726 144 : sc.push_back( "CZ" );
1727 144 : sc.push_back( "OH" );
1728 144 : sc.push_back( "HB1" );
1729 144 : sc.push_back( "HB2" );
1730 144 : sc.push_back( "HB3" );
1731 144 : sc.push_back( "HD1" );
1732 144 : sc.push_back( "HD2" );
1733 144 : sc.push_back( "HD3" );
1734 144 : sc.push_back( "HE1" );
1735 144 : sc.push_back( "HE2" );
1736 144 : sc.push_back( "HE3" );
1737 144 : sc.push_back( "HH" );
1738 144 : return sc;
1739 576 : } else if(s=="VAL") {
1740 576 : sc.push_back( "CB" );
1741 576 : sc.push_back( "CG1" );
1742 576 : sc.push_back( "CG2" );
1743 576 : sc.push_back( "HB" );
1744 576 : sc.push_back( "HG11" );
1745 576 : sc.push_back( "HG12" );
1746 576 : sc.push_back( "HG13" );
1747 576 : sc.push_back( "HG21" );
1748 576 : sc.push_back( "HG22" );
1749 576 : sc.push_back( "HG23" );
1750 576 : sc.push_back( "1HG1" );
1751 576 : sc.push_back( "2HG1" );
1752 576 : sc.push_back( "3HG1" );
1753 576 : sc.push_back( "1HG2" );
1754 576 : sc.push_back( "2HG2" );
1755 576 : sc.push_back( "3HG2" );
1756 576 : return sc;
1757 0 : } else plumed_merror("Sidechain atoms unknown: " + s);
1758 0 : }
1759 :
1760 47016 : bool CS2Backbone::isSP2(const std::string & resType, const std::string & atomName) {
1761 : bool sp2 = false;
1762 47016 : if (atomName == "C") return true;
1763 43848 : if (atomName == "O") return true;
1764 :
1765 40716 : if(resType == "TRP") {
1766 396 : if (atomName == "CG") sp2 = true;
1767 378 : else if (atomName == "CD1") sp2 = true;
1768 360 : else if (atomName == "CD2") sp2 = true;
1769 342 : else if (atomName == "CE2") sp2 = true;
1770 324 : else if (atomName == "CE3") sp2 = true;
1771 306 : else if (atomName == "CZ2") sp2 = true;
1772 288 : else if (atomName == "CZ3") sp2 = true;
1773 270 : else if (atomName == "CH2") sp2 = true;
1774 40320 : } else if (resType == "ASP") {
1775 1656 : if (atomName == "CG") sp2 = true;
1776 1494 : else if (atomName == "OD1") sp2 = true;
1777 1332 : else if (atomName == "OD2") sp2 = true;
1778 38664 : } else if (resType == "GLU") {
1779 2844 : if (atomName == "CD") sp2 = true;
1780 2628 : else if (atomName == "OE1") sp2 = true;
1781 2412 : else if (atomName == "OE2") sp2 = true;
1782 35820 : } else if (resType == "ARG") {
1783 2772 : if (atomName == "CZ") sp2 = true;
1784 33048 : } else if (resType == "HIS") {
1785 0 : if (atomName == "CG") sp2 = true;
1786 0 : else if (atomName == "ND1") sp2 = true;
1787 0 : else if (atomName == "CD2") sp2 = true;
1788 0 : else if (atomName == "CE1") sp2 = true;
1789 0 : else if (atomName == "NE2") sp2 = true;
1790 33048 : } else if (resType == "PHE") {
1791 5184 : if (atomName == "CG") sp2 = true;
1792 4896 : else if (atomName == "CD1") sp2 = true;
1793 4608 : else if (atomName == "CD2") sp2 = true;
1794 4320 : else if (atomName == "CE1") sp2 = true;
1795 4032 : else if (atomName == "CE2") sp2 = true;
1796 3744 : else if (atomName == "CZ") sp2 = true;
1797 27864 : } else if (resType == "TYR") {
1798 684 : if (atomName == "CG") sp2 = true;
1799 648 : else if (atomName == "CD1") sp2 = true;
1800 612 : else if (atomName == "CD2") sp2 = true;
1801 576 : else if (atomName == "CE1") sp2 = true;
1802 540 : else if (atomName == "CE2") sp2 = true;
1803 504 : else if (atomName == "CZ") sp2 = true;
1804 27180 : } else if (resType == "ASN") {
1805 1944 : if (atomName == "CG") sp2 = true;
1806 1782 : else if (atomName == "OD1") sp2 = true;
1807 25236 : } else if (resType == "GLN") {
1808 810 : if (atomName == "CD") sp2 = true;
1809 756 : else if (atomName == "OE1") sp2 = true;
1810 : }
1811 :
1812 : return sp2;
1813 : }
1814 :
1815 145728 : bool CS2Backbone::is_chi1_cx(const std::string & frg, const std::string & atm) {
1816 145728 : if(atm=="CG") return true;
1817 139788 : if((frg == "CYS")&&(atm =="SG")) return true;
1818 288792 : if(((frg == "ILE")||(frg == "VAL"))&&(atm == "CG1")) return true;
1819 145602 : if((frg == "SER")&&(atm == "OG")) return true;
1820 148878 : if((frg == "THR")&&(atm == "OG1")) return true;
1821 :
1822 : return false;
1823 : }
1824 :
1825 3599784 : void CS2Backbone::xdist_name_map(std::string & name) {
1826 7199568 : if((name == "OT1")||(name == "OC1")) name = "O";
1827 10799352 : else if ((name == "HN") || (name == "HT1") || (name == "H1")) name = "H";
1828 7139232 : else if ((name == "CG1")|| (name == "OG")||
1829 7159572 : (name == "SG") || (name == "OG1")) name = "CG";
1830 7040070 : else if ((name == "HA1")|| (name == "HA3")) name = "HA";
1831 3599784 : }
1832 :
1833 18 : void CS2Backbone::update() {
1834 : // write status file
1835 18 : if(getWstride()>0&& (getStep()%getWstride()==0 || getCPT()) ) writeStatus();
1836 18 : }
1837 :
1838 : }
1839 : }
|