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