LCOV - code coverage report
Current view: top level - isdb - CS2Backbone.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 916 967 94.7 %
Date: 2024-10-18 13:59:31 Functions: 22 23 95.7 %

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

Generated by: LCOV version 1.16