Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2021, Andrea Arsiccio
3 :
4 : This software is provided 'as-is', without any express or implied
5 : warranty. In no event will the authors be held liable for any damages
6 : arising from the use of this software.
7 :
8 : Permission is granted to anyone to use this software for any purpose,
9 : including commercial applications, and to alter it and redistribute it
10 : freely, subject to the following restrictions:
11 :
12 : 1. The origin of this software must not be misrepresented; you must not
13 : claim that you wrote the original software. If you use this software
14 : in a product, an acknowledgment in the product documentation would be
15 : appreciated but is not required.
16 : 2. Altered source versions must be plainly marked as such, and must not be
17 : misrepresented as being the original software.
18 : 3. This notice may not be removed or altered from any source distribution.
19 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
20 :
21 : #include "Sasa.h"
22 : #include "core/ActionRegister.h"
23 : #include "core/PlumedMain.h"
24 : #include "core/GenericMolInfo.h"
25 : #include "core/ActionSet.h"
26 : #include <cstdio>
27 : #include <iostream>
28 : #include <string>
29 : #include <cmath>
30 : #include <sstream>
31 : #include <algorithm>
32 : #include <iterator>
33 : #include <fstream>
34 :
35 :
36 : using namespace std;
37 :
38 : namespace PLMD {
39 : namespace sasa {
40 :
41 : //+PLUMEDOC SASAMOD_COLVAR SASA_HASEL
42 : /*
43 : Calculates the solvent accessible surface area (SASA) of a protein molecule, or other properties related to it.
44 :
45 : The atoms for which the SASA is desired should be indicated with the keyword ATOMS, and a pdb file of the protein must be provided in input with the MOLINFO keyword. The algorithm described in \cite Hasel1988 is used for the calculation. The radius of the solvent is assumed to be 0.14 nm, which is the radius of water molecules. Using the keyword NL_STRIDE it is also possible to specify the frequency with which the neighbor list for the calculation of SASA is updated (the default is every 10 steps).
46 :
47 : Different properties can be calculated and selected using the TYPE keyword:
48 :
49 : 1) the total SASA (TOTAL);
50 :
51 : 2) the free energy of transfer for the protein according to the transfer model (TRANSFER). This keyword can be used, for instance, to compute the transfer of a protein to different temperatures, as detailed in \cite Arsiccio-T-SASA-2021, or to different pressures, as detailed in \cite Arsiccio-P-SASA-2021, or to different osmolyte solutions, as detailed in \cite Arsiccio-C-SASA-2022.
52 :
53 :
54 : When the TRANSFER keyword is used, a file with the free energy of transfer values for the sidechains and backbone atoms should be provided (using the keyword DELTAGFILE). Such file should have the following format:
55 :
56 : \verbatim
57 :
58 : ----------------Sample DeltaG.dat file---------------------
59 : ALA 0.711019999999962
60 : ARG -2.24832799999996
61 : ASN -2.74838799999999
62 : ASP -2.5626376
63 : CYS 3.89864000000006
64 : GLN -1.76192
65 : GLU -2.38664400000002
66 : GLY 0
67 : HIS -3.58152799999999
68 : ILE 2.42634399999986
69 : LEU 1.77233599999988
70 : LYS -1.92576400000002
71 : MET -0.262827999999956
72 : PHE 1.62028800000007
73 : PRO -2.15598800000001
74 : SER -1.60934800000004
75 : THR -0.591559999999987
76 : TRP 1.22936000000027
77 : TYR 0.775547999999958
78 : VAL 2.12779200000011
79 : BACKBONE 1.00066920000002
80 : -----------------------------------------------------------
81 : \endverbatim
82 :
83 : where the second column is the free energy of transfer for each sidechain/backbone, in kJ/mol.
84 :
85 : A Python script for the computation of free energy of transfer values to describe the effect of osmolyte concentration, temperature and pressure (according to \cite Arsiccio-C-SASA-2022, \cite Arsiccio-T-SASA-2021 and \cite Arsiccio-P-SASA-2021) is freely available at https://github.com/andrea-arsiccio/DeltaG-calculation. The script automatically outputs a DeltaG.dat file compatible with this SASA module.
86 :
87 :
88 : If the DELTAGFILE is not provided, the program computes the free energy of transfer values as if they had to take into account the effect of temperature according to approaches 2 or 3 in the paper \cite Arsiccio-T-SASA-2021. Please read and cite this paper if using the transfer model for computing the effect of temperature in implicit solvent simulations. For this purpose, the keyword APPROACH should be added, and set to either 2 or 3.
89 :
90 : The SASA usually makes sense when atoms used for the calculation are all part of the same molecule. When running with periodic boundary conditions, the atoms should be in the proper periodic image. This is done automatically since PLUMED 2.2, by considering the ordered list of atoms and rebuilding the broken entities using a procedure that is equivalent to that done in \ref WHOLEMOLECULES. Notice that rebuilding is local to this action. This is different from \ref WHOLEMOLECULES which actually modifies the coordinates stored in PLUMED.
91 :
92 : In case you want to recover the old behavior you should use the NOPBC flag.
93 : In that case you need to take care that atoms are in the correct periodic image.
94 :
95 : The SASA may also be computed using the SASA_LCPO collective variable, which makes use of the LCPO algorithm \cite Weiser1999. SASA_LCPO is more accurate then SASA_HASEL, but the computation is slower.
96 :
97 :
98 : \par Examples
99 :
100 : The following input tells plumed to print the total SASA for atoms 10 to 20 in a protein chain.
101 : \plumedfile
102 : SASA_HASEL TYPE=TOTAL ATOMS=10-20 NL_STRIDE=10 LABEL=sasa
103 : PRINT ARG=sasa STRIDE=1 FILE=colvar
104 : \endplumedfile
105 :
106 :
107 : The following input tells plumed to compute the transfer free energy for the protein chain containing atoms 10 to 20. Such transfer free energy is then used as a bias in the simulation (e.g., implicit solvent simulations). The free energy of transfer values are read from a file called DeltaG.dat.
108 :
109 : \plumedfile
110 : SASA_HASEL TYPE=TRANSFER ATOMS=10-20 NL_STRIDE=10 DELTAGFILE=DeltaG.dat LABEL=sasa
111 :
112 : bias: BIASVALUE ARG=sasa
113 :
114 : PRINT ARG=sasa,bias.* STRIDE=1 FILE=colvar
115 : \endplumedfile
116 :
117 : The following input tells plumed to compute the transfer free energy for the protein chain containing atoms 10 to 20. Such transfer free energy is then used as a bias in the simulation (e.g., implicit solvent simulations). The free energy of transfer values are computed according to \cite Arsiccio-T-SASA-2021, and take into account the effect of temperature using approach 2 as described in the paper.
118 :
119 : \plumedfile
120 : SASA_HASEL TYPE=TRANSFER ATOMS=10-20 NL_STRIDE=10 APPROACH=2 LABEL=sasa
121 :
122 : bias: BIASVALUE ARG=sasa
123 :
124 : PRINT ARG=sasa,bias.* STRIDE=1 FILE=colvar
125 : \endplumedfile
126 :
127 : */
128 : //+ENDPLUMEDOC
129 :
130 : class SASA_HASEL : public Colvar {
131 : private:
132 : enum CV_TYPE {TOTAL, TRANSFER};
133 : int sasa_type;
134 : bool nopbc;
135 : double rs;
136 : string DeltaGValues;
137 : int approach;
138 : unsigned stride;
139 : unsigned nl_update;
140 : int firstStepFlag;
141 : double Ti;
142 : // cppcheck-suppress duplInheritedMember
143 : std::vector<AtomNumber> atoms;
144 : vector < vector < std::string > > AtomResidueName;
145 : vector < vector < double > > SASAparam;
146 : vector < vector < std::string > > CONNECTparam;
147 : unsigned natoms;
148 : vector < vector < double > > MaxSurf;
149 : vector < vector < double > > DeltaG;
150 : vector < vector < int > > Nlist;
151 : public:
152 : static void registerKeywords(Keywords& keys);
153 : explicit SASA_HASEL(const ActionOptions&);
154 : void readPDB();
155 : map<string, vector<std::string> > setupHASELparam();
156 : void readSASAparam();
157 : void calcNlist();
158 : map<string, vector<double> > setupMaxSurfMap();
159 : void readMaxSurf();
160 : void readDeltaG();
161 : void computeDeltaG();
162 : virtual void calculate();
163 : };
164 :
165 : PLUMED_REGISTER_ACTION(SASA_HASEL,"SASA_HASEL")
166 :
167 4 : void SASA_HASEL::registerKeywords(Keywords& keys) {
168 4 : Colvar::registerKeywords(keys);
169 8 : keys.add("atoms","ATOMS","the group of atoms that you are calculating the SASA for");
170 8 : keys.add("compulsory","TYPE","TOTAL","The type of calculation you want to perform. Can be TOTAL or TRANSFER");
171 8 : keys.add("compulsory", "NL_STRIDE", "The frequency with which the neighbor list for the calculation of SASA is updated.");
172 8 : keys.add("optional","DELTAGFILE","a file containing the free energy of transfer values for backbone and sidechains atoms. Necessary only if TYPE = TRANSFER. A Python script for the computation of free energy of transfer values to describe the effect of osmolyte concentration, temperature and pressure is freely available at https://github.com/andrea-arsiccio/DeltaG-calculation. The script automatically outputs a DeltaG.dat file compatible with this SASA module. If TYPE = TRANSFER and no DELTAGFILE is provided, the free energy values are those describing the effect of temperature, and are computed using the temperature value passed by the MD engine");
173 8 : keys.add("optional","APPROACH","either approach 2 or 3. Necessary only if TYPE = TRANSFER and no DELTAGFILE is provided. If TYPE = TRANSFER and no DELTAGFILE is provided, the free energy values are those describing the effect of temperature, and the program must know if approach 2 or 3 (as described in Arsiccio and Shea, Protein Cold Denaturation in Implicit Solvent Simulations: A Transfer Free Energy Approach, J. Phys. Chem. B, 2021) needs to be used to compute them");
174 8 : keys.setValueDescription("scalar","the solvent accessible surface area (SASA) of the molecule");
175 4 : }
176 :
177 :
178 2 : SASA_HASEL::SASA_HASEL(const ActionOptions&ao):
179 : PLUMED_COLVAR_INIT(ao),
180 2 : nopbc(false),
181 2 : stride(10),
182 2 : nl_update(0),
183 4 : DeltaGValues("absent"),
184 2 : Ti(0),
185 2 : firstStepFlag(0)
186 : {
187 2 : rs = 0.14;
188 2 : parse("DELTAGFILE",DeltaGValues);
189 2 : parse("APPROACH", approach);
190 4 : parseAtomList("ATOMS",atoms);
191 2 : if(atoms.size()==0) error("no atoms specified");
192 : std::string Type;
193 2 : parse("TYPE",Type);
194 2 : parse("NL_STRIDE", stride);
195 2 : parseFlag("NOPBC",nopbc);
196 2 : checkRead();
197 :
198 2 : if(Type=="TOTAL") sasa_type=TOTAL;
199 1 : else if(Type=="TRANSFER") sasa_type=TRANSFER;
200 0 : else error("Unknown SASA type");
201 :
202 2 : switch(sasa_type)
203 : {
204 1 : case TOTAL: log.printf(" TOTAL SASA;"); break;
205 1 : case TRANSFER: log.printf(" TRANSFER MODEL;"); break;
206 : }
207 :
208 2 : log.printf(" atoms involved : ");
209 242 : for(unsigned i=0; i<atoms.size(); ++i) {
210 240 : if(i%25==0) log<<"\n";
211 240 : log.printf("%d ",atoms[i].serial());
212 : }
213 2 : log.printf("\n");
214 :
215 2 : if(nopbc) {
216 0 : log<<" PBC will be ignored\n";
217 : } else {
218 2 : log<<" broken molecules will be rebuilt assuming atoms are in the proper order\n";
219 : }
220 :
221 :
222 4 : addValueWithDerivatives(); setNotPeriodic();
223 2 : requestAtoms(atoms);
224 :
225 2 : natoms = getNumberOfAtoms();
226 2 : AtomResidueName.resize(2);
227 2 : SASAparam.resize(natoms);
228 2 : CONNECTparam.resize(natoms);
229 2 : MaxSurf.resize(natoms);
230 2 : DeltaG.resize(natoms+1);
231 2 : Nlist.resize(natoms);
232 :
233 :
234 2 : }
235 :
236 :
237 : //splits strings into tokens. Used to read into SASA parameters file and into reference pdb file
238 : template <class Container>
239 42 : void split(const std::string& str, Container& cont)
240 : {
241 42 : std::istringstream iss(str);
242 84 : std::copy(std::istream_iterator<std::string>(iss),
243 42 : std::istream_iterator<std::string>(),
244 : std::back_inserter(cont));
245 42 : }
246 :
247 :
248 : //reads input PDB file provided with MOLINFO, and assigns atom and residue names to each atom involved in the CV
249 :
250 2 : void SASA_HASEL::readPDB() {
251 2 : auto* moldat = plumed.getActionSet().selectLatest<GenericMolInfo*>(this);
252 2 : if( ! moldat ) error("Unable to find MOLINFO in input");
253 : AtomResidueName[0].clear();
254 : AtomResidueName[1].clear();
255 :
256 242 : for(unsigned i=0; i<natoms; i++) {
257 240 : string Aname = moldat->getAtomName(atoms[i]);
258 240 : string Rname = moldat->getResidueName(atoms[i]);
259 240 : AtomResidueName[0].push_back (Aname);
260 240 : AtomResidueName[1].push_back (Rname);
261 : }
262 :
263 2 : }
264 :
265 :
266 : //Hasel et al. parameters database
267 2 : map<string, vector<std::string> > SASA_HASEL::setupHASELparam() {
268 : map<string, vector<std::string> > haselmap;
269 16 : haselmap["ALA_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
270 16 : haselmap["ALA_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
271 16 : haselmap["ALA_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
272 16 : haselmap["ALA_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
273 16 : haselmap["ALA_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
274 16 : haselmap["ALA_CB"] = { "2.0", "0.88", "CA", "Z", "Z", "Z", };
275 16 : haselmap["ASP_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
276 16 : haselmap["ASP_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
277 16 : haselmap["ASP_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
278 16 : haselmap["ASP_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
279 16 : haselmap["ASP_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
280 16 : haselmap["ASP_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
281 16 : haselmap["ASP_CG"] = { "1.72", "1.554", "CB", "OD1", "OD2", "Z", };
282 16 : haselmap["ASP_OD1"] = { "1.5", "0.926", "CG", "Z", "Z", "Z", };
283 16 : haselmap["ASP_OD2"] = { "1.7", "0.922", "CG", "Z", "Z", "Z", };
284 16 : haselmap["ASN_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
285 16 : haselmap["ASN_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
286 16 : haselmap["ASN_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
287 16 : haselmap["ASN_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
288 16 : haselmap["ASN_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
289 16 : haselmap["ASN_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
290 16 : haselmap["ASN_CG"] = { "1.7", "2.149", "CB", "OD1", "ND2", "Z", };
291 16 : haselmap["ASN_OD1"] = { "1.5", "0.926", "CG", "Z", "Z", "Z", };
292 16 : haselmap["ASN_ND2"] = { "1.6", "1.215", "CG", "1HD2", "1HD2", "Z", };
293 16 : haselmap["ASN_1HD2"] = { "1.1", "1.128", "ND2", "Z", "Z", "Z", };
294 16 : haselmap["ASN_2HD2"] = { "1.1", "1.128", "ND2", "Z", "Z", "Z", };
295 16 : haselmap["ARG_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
296 16 : haselmap["ARG_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
297 16 : haselmap["ARG_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
298 16 : haselmap["ARG_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
299 16 : haselmap["ARG_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
300 16 : haselmap["ARG_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
301 16 : haselmap["ARG_CG"] = { "1.9", "1.045", "CB", "CD", "Z", "Z", };
302 16 : haselmap["ARG_CD"] = { "1.9", "1.045", "CG", "NE", "Z", "Z", };
303 16 : haselmap["ARG_NE"] = { "1.55", "1.028", "CD", "HE", "CZ", "Z", };
304 16 : haselmap["ARG_NH1"] = { "1.55", "1.028", "CZ", "1HH1", "2HH1", "Z", };
305 16 : haselmap["ARG_NH2"] = { "1.55", "1.028", "CZ", "1HH2", "2HH2", "Z", };
306 16 : haselmap["ARG_CZ"] = { "1.72", "1.554", "NE", "NH1", "NH2", "Z", };
307 16 : haselmap["ARG_HE"] = { "1.1", "1.128", "NE", "Z", "Z", "Z", };
308 16 : haselmap["ARG_1HH2"] = { "1.1", "1.128", "NH2", "Z", "Z", "Z", };
309 16 : haselmap["ARG_2HH2"] = { "1.1", "1.128", "NH2", "Z", "Z", "Z", };
310 16 : haselmap["ARG_2HH1"] = { "1.1", "1.128", "NH1", "Z", "Z", "Z", };
311 16 : haselmap["ARG_1HH1"] = { "1.1", "1.128", "NH1", "Z", "Z", "Z", };
312 16 : haselmap["CYS_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
313 16 : haselmap["CYS_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
314 16 : haselmap["CYS_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
315 16 : haselmap["CYS_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
316 16 : haselmap["CYS_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
317 16 : haselmap["CYS_CB"] = { "1.9", "1.045", "CA", "SG", "Z", "Z", };
318 16 : haselmap["CYS_SG"] = { "1.8", "1.121", "CB", "HG", "Z", "Z", };
319 16 : haselmap["CYS_HG"] = { "1.2", "0.928", "SG", "Z", "Z", "Z", };
320 16 : haselmap["GLU_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
321 16 : haselmap["GLU_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
322 16 : haselmap["GLU_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
323 16 : haselmap["GLU_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
324 16 : haselmap["GLU_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
325 16 : haselmap["GLU_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
326 16 : haselmap["GLU_CG"] = { "1.9", "1.045", "CB", "CD", "Z", "Z", };
327 16 : haselmap["GLU_CD"] = { "1.72", "1.554", "CG", "OE1", "OE2", "Z", };
328 16 : haselmap["GLU_OE1"] = { "1.5", "0.926", "CD", "Z", "Z", "Z", };
329 16 : haselmap["GLU_OE2"] = { "1.7", "0.922", "CD", "Z", "Z", "Z", };
330 16 : haselmap["GLN_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
331 16 : haselmap["GLN_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
332 16 : haselmap["GLN_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
333 16 : haselmap["GLN_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
334 16 : haselmap["GLN_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
335 16 : haselmap["GLN_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
336 16 : haselmap["GLN_CG"] = { "1.9", "1.045", "CB", "CD", "Z", "Z", };
337 16 : haselmap["GLN_CD"] = { "1.72", "1.554", "CG", "OE1", "NE2", "Z", };
338 16 : haselmap["GLN_OE1"] = { "1.5", "0.926", "CD", "Z", "Z", "Z", };
339 16 : haselmap["GLN_NE2"] = { "1.6", "1.215", "CD", "2HE2", "1HE2", "Z", };
340 16 : haselmap["GLN_2HE2"] = { "1.1", "1.128", "NE2", "Z", "Z", "Z", };
341 16 : haselmap["GLN_1HE2"] = { "1.1", "1.128", "NE2", "Z", "Z", "Z", };
342 16 : haselmap["GLY_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
343 16 : haselmap["GLY_CA"] = { "1.7", "2.149", "N", "C", "Z", "Z", };
344 16 : haselmap["GLY_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
345 16 : haselmap["GLY_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
346 16 : haselmap["GLY_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
347 16 : haselmap["HIS_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
348 16 : haselmap["HIS_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
349 16 : haselmap["HIS_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
350 16 : haselmap["HIS_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
351 16 : haselmap["HIS_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
352 16 : haselmap["HIS_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
353 16 : haselmap["HIS_CG"] = { "1.72", "1.554", "CB", "CD2", "ND1", "Z", };
354 16 : haselmap["HIS_ND1"] = { "1.55", "1.028", "CG", "CE1", "Z", "Z", };
355 16 : haselmap["HIS_CE1"] = { "1.8", "1.073", "ND1", "NE2", "Z", "Z", };
356 16 : haselmap["HIS_NE2"] = { "1.55", "1.413", "CE1", "2HE", "CD2", "Z", };
357 16 : haselmap["HIS_CD2"] = { "1.8", "1.073", "NE2", "CG", "Z", "Z", };
358 16 : haselmap["HIS_2HE"] = { "1.1", "1.128", "NE2", "Z", "Z", "Z", };
359 16 : haselmap["ILE_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
360 16 : haselmap["ILE_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
361 16 : haselmap["ILE_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
362 16 : haselmap["ILE_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
363 16 : haselmap["ILE_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
364 16 : haselmap["ILE_CB"] = { "1.8", "1.276", "CA", "CG2", "CG1", "Z", };
365 16 : haselmap["ILE_CG2"] = { "2.0", "0.88", "CB", "Z", "Z", "Z", };
366 16 : haselmap["ILE_CG1"] = { "1.9", "1.045", "CB", "CD1", "Z", "Z", };
367 16 : haselmap["ILE_CD1"] = { "2.0", "0.88", "CG1", "Z", "Z", "Z", };
368 16 : haselmap["LEU_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
369 16 : haselmap["LEU_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
370 16 : haselmap["LEU_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
371 16 : haselmap["LEU_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
372 16 : haselmap["LEU_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
373 16 : haselmap["LEU_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
374 16 : haselmap["LEU_CG"] = { "1.8", "1.276", "CB", "CD1", "CD2", "Z", };
375 16 : haselmap["LEU_CD1"] = { "2.0", "0.88", "CG", "Z", "Z", "Z", };
376 16 : haselmap["LEU_CD2"] = { "2.0", "0.88", "CG", "Z", "Z", "Z", };
377 16 : haselmap["LYS_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
378 16 : haselmap["LYS_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
379 16 : haselmap["LYS_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
380 16 : haselmap["LYS_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
381 16 : haselmap["LYS_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
382 16 : haselmap["LYS_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
383 16 : haselmap["LYS_CG"] = { "1.9", "1.045", "CB", "CD", "Z", "Z", };
384 16 : haselmap["LYS_CD"] = { "1.9", "1.045", "CG", "CE", "Z", "Z", };
385 16 : haselmap["LYS_CE"] = { "1.9", "1.045", "CD", "NZ", "Z", "Z", };
386 16 : haselmap["LYS_NZ"] = { "1.6", "1.215", "CE", "1HZ", "2HZ", "3HZ", };
387 16 : haselmap["LYS_1HZ"] = { "1.1", "1.128", "NZ", "Z", "Z", "Z", };
388 16 : haselmap["LYS_2HZ"] = { "1.1", "1.128", "NZ", "Z", "Z", "Z", };
389 16 : haselmap["LYS_3HZ"] = { "1.1", "1.128", "NZ", "Z", "Z", "Z", };
390 16 : haselmap["MET_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
391 16 : haselmap["MET_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
392 16 : haselmap["MET_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
393 16 : haselmap["MET_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
394 16 : haselmap["MET_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
395 16 : haselmap["MET_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
396 16 : haselmap["MET_CG"] = { "1.9", "1.045", "CB", "SD", "Z", "Z", };
397 16 : haselmap["MET_SD"] = { "1.8", "1.121", "CG", "CE", "Z", "Z", };
398 16 : haselmap["MET_CE"] = { "2.0", "0.88", "SD", "Z", "Z", "Z", };
399 16 : haselmap["PHE_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
400 16 : haselmap["PHE_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
401 16 : haselmap["PHE_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
402 16 : haselmap["PHE_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
403 16 : haselmap["PHE_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
404 16 : haselmap["PHE_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
405 16 : haselmap["PHE_CG"] = { "1.72", "1.554", "CB", "CD1", "CD2", "Z", };
406 16 : haselmap["PHE_CD1"] = { "1.8", "1.073", "CG", "CE1", "Z", "Z", };
407 16 : haselmap["PHE_CE1"] = { "1.8", "1.073", "CD1", "CZ", "Z", "Z", };
408 16 : haselmap["PHE_CZ"] = { "1.8", "1.073", "CE1", "CE2", "Z", "Z", };
409 16 : haselmap["PHE_CE2"] = { "1.8", "1.073", "CZ", "CD2", "Z", "Z", };
410 16 : haselmap["PHE_CD2"] = { "1.8", "1.073", "CE2", "CG", "Z", "Z", };
411 16 : haselmap["PRO_N"] = { "1.55", "1.028", "CD", "CA", "Z", "Z", };
412 16 : haselmap["PRO_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
413 16 : haselmap["PRO_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
414 16 : haselmap["PRO_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
415 16 : haselmap["PRO_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
416 16 : haselmap["PRO_CG"] = { "1.9", "1.045", "CB", "CD", "Z", "Z", };
417 16 : haselmap["PRO_CD"] = { "1.9", "1.045", "CG", "N", "Z", "Z", };
418 16 : haselmap["SER_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
419 16 : haselmap["SER_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
420 16 : haselmap["SER_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
421 16 : haselmap["SER_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
422 16 : haselmap["SER_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
423 16 : haselmap["SER_CB"] = { "1.9", "1.045", "CA", "OG", "Z", "Z", };
424 16 : haselmap["SER_OG"] = { "1.52", "1.08", "CB", "HG", "Z", "Z", };
425 16 : haselmap["SER_HG"] = { "1.0", "0.944", "OG", "Z", "Z", "Z", };
426 16 : haselmap["THR_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
427 16 : haselmap["THR_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
428 16 : haselmap["THR_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
429 16 : haselmap["THR_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
430 16 : haselmap["THR_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
431 16 : haselmap["THR_CB"] = { "1.8", "1.276", "CA", "CG2", "OG1", "Z", };
432 16 : haselmap["THR_CG2"] = { "2.0", "0.88", "CB", "Z", "Z", "Z", };
433 16 : haselmap["THR_OG1"] = { "1.52", "1.08", "1HG", "CB", "Z", "Z", };
434 16 : haselmap["THR_1HG"] = { "1.0", "0.944", "OG1", "Z", "Z", "Z", };
435 16 : haselmap["TRP_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
436 16 : haselmap["TRP_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
437 16 : haselmap["TRP_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
438 16 : haselmap["TRP_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
439 16 : haselmap["TRP_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
440 16 : haselmap["TRP_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
441 16 : haselmap["TRP_CG"] = { "1.72", "1.554", "CB", "CD2", "CD1", "Z", };
442 16 : haselmap["TRP_CD1"] = { "1.8", "1.073", "CG", "NE1", "Z", "Z", };
443 16 : haselmap["TRP_NE1"] = { "1.55", "1.413", "CD1", "CE2", "1HE", "Z", };
444 16 : haselmap["TRP_CE2"] = { "1.72", "1.554", "NE1", "CD2", "CZ2", "Z", };
445 16 : haselmap["TRP_CZ2"] = { "1.8", "1.073", "CE2", "CH2", "Z", "Z", };
446 16 : haselmap["TRP_CH2"] = { "1.8", "1.073", "CZ2", "CZ3", "Z", "Z", };
447 16 : haselmap["TRP_CZ3"] = { "1.8", "1.073", "CH2", "CE3", "Z", "Z", };
448 16 : haselmap["TRP_CE3"] = { "1.8", "1.073", "CZ3", "CD2", "Z", "Z", };
449 16 : haselmap["TRP_CD2"] = { "1.72", "1.554", "CE3", "CE2", "CG", "Z", };
450 16 : haselmap["TRP_1HE"] = { "1.1", "1.128", "NE1", "Z", "Z", "Z", };
451 16 : haselmap["TYR_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
452 16 : haselmap["TYR_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
453 16 : haselmap["TYR_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
454 16 : haselmap["TYR_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
455 16 : haselmap["TYR_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
456 16 : haselmap["TYR_CB"] = { "1.9", "1.045", "CA", "CG", "Z", "Z", };
457 16 : haselmap["TYR_CG"] = { "1.72", "1.554", "CB", "CD1", "CD2", "Z", };
458 16 : haselmap["TYR_CD1"] = { "1.8", "1.073", "CG", "CE1", "Z", "Z", };
459 16 : haselmap["TYR_CE1"] = { "1.8", "1.073", "CD1", "CZ", "Z", "Z", };
460 16 : haselmap["TYR_CZ"] = { "1.72", "1.554", "CE1", "OH", "CE2", "Z", };
461 16 : haselmap["TYR_OH"] = { "1.52", "1.08", "CZ", "HH", "Z", "Z", };
462 16 : haselmap["TYR_HH"] = { "1.0", "0.944", "OH", "Z", "Z", "Z", };
463 16 : haselmap["TYR_CE2"] = { "1.8", "1.073", "CZ", "CD2", "Z", "Z", };
464 16 : haselmap["TYR_CD2"] = { "1.8", "1.073", "CE2", "CG", "Z", "Z", };
465 16 : haselmap["VAL_N"] = { "1.55", "1.028", "H", "CA", "Z", "Z", };
466 16 : haselmap["VAL_CA"] = { "1.7", "2.149", "N", "C", "CB", "Z", };
467 16 : haselmap["VAL_C"] = { "1.72", "1.554", "CA", "O", "Z", "Z", };
468 16 : haselmap["VAL_O"] = { "1.5", "0.926", "C", "Z", "Z", "Z", };
469 16 : haselmap["VAL_H"] = { "1.1", "1.128", "N", "Z", "Z", "Z", };
470 16 : haselmap["VAL_CB"] = { "1.8", "1.276", "CA", "CG1", "CG2", "Z", };
471 16 : haselmap["VAL_CG1"] = { "2.0", "0.88", "CB", "Z", "Z", "Z", };
472 16 : haselmap["VAL_CG2"] = { "2.0", "0.88", "CB", "Z", "Z", "Z", };
473 2 : return haselmap;
474 : }
475 :
476 : //assigns SASA parameters to each atom reading from HASEL parameter database
477 2 : void SASA_HASEL::readSASAparam() {
478 :
479 242 : for(unsigned i=0; i<natoms; i++) {
480 240 : SASAparam[i].clear();
481 : CONNECTparam[i].clear();
482 : }
483 :
484 : map<string, vector<std::string> > haselmap;
485 4 : haselmap = setupHASELparam();
486 : vector<std::string> HASELparamVector;
487 : string identifier;
488 :
489 :
490 242 : for(unsigned i=0; i<natoms; i++) {
491 480 : identifier = AtomResidueName[1][i]+"_"+AtomResidueName[0][i];
492 240 : if (haselmap.find(identifier)!=haselmap.end()) {
493 144 : HASELparamVector = haselmap.at(identifier);
494 144 : SASAparam[i].push_back (std::atof(HASELparamVector[0].c_str())+rs*10);
495 144 : SASAparam[i].push_back (std::atof(HASELparamVector[1].c_str()));
496 288 : CONNECTparam[i].push_back (HASELparamVector[2].c_str());
497 288 : CONNECTparam[i].push_back (HASELparamVector[3].c_str());
498 288 : CONNECTparam[i].push_back (HASELparamVector[4].c_str());
499 288 : CONNECTparam[i].push_back (HASELparamVector[5].c_str());
500 : }
501 : }
502 :
503 :
504 242 : for(unsigned i=0; i<natoms; i++) {
505 240 : if (SASAparam[i].size()==0 ) {
506 96 : if ((AtomResidueName[0][i][0]=='O') || (AtomResidueName[0][i][0]=='N') || (AtomResidueName[0][i][0]=='C') || (AtomResidueName[0][i][0]=='S')) {
507 0 : cout << "Could not find SASA paramaters for atom " << AtomResidueName[0][i] << " of residue " << AtomResidueName[1][i] << endl;
508 0 : error ("missing SASA parameters\n");
509 : }
510 : }
511 : }
512 :
513 :
514 4 : }
515 :
516 :
517 :
518 : //Max Surf values, used only if TYPE=TRANSFER
519 1 : map<string, vector<double> > SASA_HASEL::setupMaxSurfMap() {
520 : // Max Surface Area for backbone and sidechain, in nm2
521 : map<string, vector<double> > valuemap;
522 1 : valuemap ["ALA"]= {0.56425,0.584851,};
523 1 : valuemap ["ARG"]= {0.498656,1.808093,};
524 1 : valuemap ["ASN"]= {0.473409,0.818394,};
525 1 : valuemap ["ASP"]= {0.477057,0.977303,};
526 1 : valuemap ["CYS"]= {0.507512,0.791483,};
527 1 : valuemap ["GLN"]= {0.485859,1.281534,};
528 1 : valuemap ["GLU"]= {0.495054,1.464718,};
529 1 : valuemap ["GLY"]= {0.658632,0,};
530 1 : valuemap ["HIS"]= {0.48194,1.118851,};
531 1 : valuemap ["ILE"]= {0.461283,1.450569,};
532 1 : valuemap ["LEU"]= {0.476315,1.498843,};
533 1 : valuemap ["LYS"]= {0.493533,1.619731,};
534 1 : valuemap ["MET"]= {0.507019,1.631904,};
535 1 : valuemap ["PHE"]= {0.457462, 1.275125,};
536 1 : valuemap ["PRO"]= {0.315865,0.859456,};
537 1 : valuemap ["SER"]= {0.48636,0.627233,};
538 1 : valuemap ["THR"]= {0.45064,0.91088,};
539 1 : valuemap ["TRP"]= {0.45762,1.366369,};
540 1 : valuemap ["TYR"]= {0.461826,1.425822,};
541 1 : valuemap ["VAL"]= {0.477054,1.149101,};
542 1 : return valuemap;
543 : }
544 :
545 :
546 :
547 : //reads maximum surface values per residue type and assigns values to each atom, only used if sasa_type = TRANSFER
548 :
549 1 : void SASA_HASEL::readMaxSurf() {
550 : map<string, vector<double> > valuemap;
551 2 : valuemap = setupMaxSurfMap();
552 : vector<double> MaxSurfVector;
553 :
554 121 : for(unsigned i=0; i<natoms; i++) {
555 120 : MaxSurf[i].clear();
556 120 : MaxSurfVector = valuemap.at(AtomResidueName[1][i]);
557 120 : MaxSurf[i].push_back (MaxSurfVector[0]*100);
558 120 : MaxSurf[i].push_back (MaxSurfVector[1]*100);
559 : }
560 1 : }
561 :
562 : //reads file with free energy values for sidechains and for the backbone, and assigns values to each atom. Only used if sasa_type = TRANSFER
563 :
564 1 : void SASA_HASEL::readDeltaG() {
565 :
566 121 : for(unsigned i=0; i<natoms; i++) {
567 120 : DeltaG[i].clear();
568 : }
569 :
570 : string DeltaGline;
571 1 : fstream DeltaGFile;
572 1 : DeltaGFile.open(DeltaGValues);
573 1 : if (DeltaGFile) {
574 : int backboneflag = 0;
575 23 : while(getline(DeltaGFile, DeltaGline)) {
576 22 : if(!DeltaGline.empty()) {
577 : std::vector<std::string> DeltaGtoken;
578 21 : split(DeltaGline, DeltaGtoken);
579 2541 : for(unsigned i=0; i<natoms; i++) {
580 2520 : if (DeltaGtoken[0].compare(AtomResidueName[1][i])==0 ) {
581 120 : DeltaG[i].push_back (std::atof(DeltaGtoken[1].c_str()));
582 : }
583 : }
584 21 : if (DeltaGtoken[0].compare("BACKBONE")==0 ) {
585 : backboneflag = 1;
586 1 : DeltaG[natoms].push_back (std::atof(DeltaGtoken[1].c_str()));
587 : }
588 21 : }
589 : }
590 1 : if ( backboneflag == 0) error("Cannot find backbone value in Delta G parameters file\n");
591 : }
592 0 : else error("Cannot open DeltaG file");
593 :
594 121 : for(unsigned i=0; i<natoms; i++) {
595 120 : if (DeltaG[i].size()==0 ) {
596 0 : cout << "Delta G value for residue " << AtomResidueName[1][i] << " not found " << endl;
597 0 : error ("missing Delta G parameter\n");
598 : }
599 : }
600 :
601 2 : }
602 :
603 : //computes free energy values for the sidechains and for the backbone, and assigns values to each atom. Only used if sasa_type = TRANSFER, and if no DELTAGFILE is provided. In this case, the free energy values are those describing the effect of temperature, and the program must know if approach 2 or 3 (as described in Arsiccio and Shea, Protein Cold Denaturation in Implicit Solvent Simulations: A Transfer Free Energy Approach, JPCB, 2021) needs to be used to compute them.
604 :
605 0 : void SASA_HASEL::computeDeltaG() {
606 :
607 0 : for(unsigned i=0; i<natoms; i++) {
608 0 : DeltaG[i].clear();
609 : }
610 :
611 : double T;
612 0 : T = getkBT()/getKBoltzmann();
613 :
614 0 : if (T != Ti) {
615 0 : for(unsigned i=0; i<natoms; i++) {
616 0 : if (approach==2) {
617 0 : if (AtomResidueName[1][i].compare("ALA")==0) {
618 0 : DeltaG[i].push_back (-2.995/1000*std::pow(T,2)+1.808*T-272.895);
619 : }
620 0 : if (AtomResidueName[1][i].compare("ARG")==0) {
621 0 : DeltaG[i].push_back (-3.182/1000*std::pow(T,2)+1.894*T-282.032);
622 : }
623 0 : if (AtomResidueName[1][i].compare("ASN")==0) {
624 0 : DeltaG[i].push_back (-1.047/1000*std::pow(T,2)+0.6068*T-87.846);
625 : }
626 0 : if (AtomResidueName[1][i].compare("ASP")==0) {
627 0 : DeltaG[i].push_back (-0.1794/1000*std::pow(T,2)+0.1091*T-16.526);
628 : }
629 0 : if (AtomResidueName[1][i].compare("CYS")==0) {
630 0 : DeltaG[i].push_back (-3.09/1000*std::pow(T,2)+1.835*T-272.26);
631 : }
632 0 : if (AtomResidueName[1][i].compare("GLN")==0) {
633 0 : DeltaG[i].push_back (-2.23/1000*std::pow(T,2)+1.335*T-199.707);
634 : }
635 0 : if (AtomResidueName[1][i].compare("GLU")==0) {
636 0 : DeltaG[i].push_back (-1.511/1000*std::pow(T,2)+0.8904*T-131.168);
637 : }
638 0 : if (AtomResidueName[1][i].compare("GLY")==0) {
639 0 : DeltaG[i].push_back (0);
640 : }
641 0 : if (AtomResidueName[1][i].compare("HIS")==0) {
642 0 : DeltaG[i].push_back (-3.482/1000*std::pow(T,2)+2.084*T-311.694);
643 : }
644 0 : if (AtomResidueName[1][i].compare("ILE")==0) {
645 0 : DeltaG[i].push_back (-6.364/1000*std::pow(T,2)+3.8*T-567.444);
646 : }
647 0 : if (AtomResidueName[1][i].compare("LEU")==0) {
648 0 : DeltaG[i].push_back (-7.466/1000*std::pow(T,2)+4.417*T-653.394);
649 : }
650 0 : if (AtomResidueName[1][i].compare("LYS")==0) {
651 0 : DeltaG[i].push_back (-2.091/1000*std::pow(T,2)+1.2458*T-185.549);
652 : }
653 0 : if (AtomResidueName[1][i].compare("MET")==0) {
654 0 : DeltaG[i].push_back (-3.807/1000*std::pow(T,2)+2.272*T-339.007);
655 : }
656 0 : if (AtomResidueName[1][i].compare("PHE")==0) {
657 0 : DeltaG[i].push_back (-7.828/1000*std::pow(T,2)+4.644*T-688.874);
658 : }
659 0 : if (AtomResidueName[1][i].compare("PRO")==0) {
660 0 : DeltaG[i].push_back (-2.347/1000*std::pow(T,2)+1.411*T-212.059);
661 : }
662 0 : if (AtomResidueName[1][i].compare("SER")==0) {
663 0 : DeltaG[i].push_back (1.813/1000*std::pow(T,2)-1.05*T+151.957);
664 : }
665 0 : if (AtomResidueName[1][i].compare("THR")==0) {
666 0 : DeltaG[i].push_back (-2.64/1000*std::pow(T,2)+1.591*T-239.516);
667 : }
668 0 : if (AtomResidueName[1][i].compare("TRP")==0) {
669 0 : DeltaG[i].push_back (-11.66/1000*std::pow(T,2)+6.916*T-1025.293);
670 : }
671 0 : if (AtomResidueName[1][i].compare("TYR")==0) {
672 0 : DeltaG[i].push_back (-7.513/1000*std::pow(T,2)+4.478*T-667.261);
673 : }
674 0 : if (AtomResidueName[1][i].compare("VAL")==0) {
675 0 : DeltaG[i].push_back (-4.902/1000*std::pow(T,2)+2.921*T-435.309);
676 : }
677 0 : DeltaG[natoms].push_back (-0.6962/1000*std::pow(T,2)+0.4426*T-70.015);
678 : }
679 0 : if (approach==3) {
680 0 : if (AtomResidueName[1][i].compare("ALA")==0) {
681 0 : DeltaG[i].push_back (-2.995/1000*std::pow(T,2)+1.808*T-272.105);
682 : }
683 0 : if (AtomResidueName[1][i].compare("ARG")==0) {
684 0 : DeltaG[i].push_back (-3.182/1000*std::pow(T,2)+1.894*T-284.086);
685 : }
686 0 : if (AtomResidueName[1][i].compare("ASN")==0) {
687 0 : DeltaG[i].push_back (-1.047/1000*std::pow(T,2)+0.6068*T-90.597);
688 : }
689 0 : if (AtomResidueName[1][i].compare("ASP")==0) {
690 0 : DeltaG[i].push_back (-0.1794/1000*std::pow(T,2)+0.1091*T-19.143);
691 : }
692 0 : if (AtomResidueName[1][i].compare("CYS")==0) {
693 0 : DeltaG[i].push_back (-3.09/1000*std::pow(T,2)+1.835*T-268.527);
694 : }
695 0 : if (AtomResidueName[1][i].compare("GLN")==0) {
696 0 : DeltaG[i].push_back (-2.23/1000*std::pow(T,2)+1.335*T-201.559);
697 : }
698 0 : if (AtomResidueName[1][i].compare("GLU")==0) {
699 0 : DeltaG[i].push_back (-1.511/1000*std::pow(T,2)+0.8904*T-133.543);
700 : }
701 0 : if (AtomResidueName[1][i].compare("GLY")==0) {
702 0 : DeltaG[i].push_back (0);
703 : }
704 0 : if (AtomResidueName[1][i].compare("HIS")==0) {
705 0 : DeltaG[i].push_back (-3.482/1000*std::pow(T,2)+2.084*T-315.398);
706 : }
707 0 : if (AtomResidueName[1][i].compare("ILE")==0) {
708 0 : DeltaG[i].push_back (-6.364/1000*std::pow(T,2)+3.8*T-564.825);
709 : }
710 0 : if (AtomResidueName[1][i].compare("LEU")==0) {
711 0 : DeltaG[i].push_back (-7.466/1000*std::pow(T,2)+4.417*T-651.483);
712 : }
713 0 : if (AtomResidueName[1][i].compare("LYS")==0) {
714 0 : DeltaG[i].push_back (-2.091/1000*std::pow(T,2)+1.2458*T-187.485);
715 : }
716 0 : if (AtomResidueName[1][i].compare("MET")==0) {
717 0 : DeltaG[i].push_back (-3.807/1000*std::pow(T,2)+2.272*T-339.242);
718 : }
719 0 : if (AtomResidueName[1][i].compare("PHE")==0) {
720 0 : DeltaG[i].push_back (-7.828/1000*std::pow(T,2)+4.644*T-687.134);
721 : }
722 0 : if (AtomResidueName[1][i].compare("PRO")==0) {
723 0 : DeltaG[i].push_back (-2.347/1000*std::pow(T,2)+1.411*T-214.211);
724 : }
725 0 : if (AtomResidueName[1][i].compare("SER")==0) {
726 0 : DeltaG[i].push_back (1.813/1000*std::pow(T,2)-1.05*T+150.289);
727 : }
728 0 : if (AtomResidueName[1][i].compare("THR")==0) {
729 0 : DeltaG[i].push_back (-2.64/1000*std::pow(T,2)+1.591*T-240.267);
730 : }
731 0 : if (AtomResidueName[1][i].compare("TRP")==0) {
732 0 : DeltaG[i].push_back (-11.66/1000*std::pow(T,2)+6.916*T-1024.284);
733 : }
734 0 : if (AtomResidueName[1][i].compare("TYR")==0) {
735 0 : DeltaG[i].push_back (-7.513/1000*std::pow(T,2)+4.478*T-666.484);
736 : }
737 0 : if (AtomResidueName[1][i].compare("VAL")==0) {
738 0 : DeltaG[i].push_back (-4.902/1000*std::pow(T,2)+2.921*T-433.013);
739 : }
740 0 : DeltaG[natoms].push_back (-0.6927/1000*std::pow(T,2)+0.4404*T-68.724);
741 : }
742 : }
743 : }
744 :
745 0 : Ti = T;
746 :
747 0 : if (firstStepFlag ==0) {
748 0 : if (approach!=2 && approach!=3) {
749 0 : cout << "Unknown approach " << approach << endl;
750 : }
751 0 : for(unsigned i=0; i<natoms; i++) {
752 0 : if (DeltaG[i].size()==0 ) {
753 0 : cout << "Delta G value for residue " << AtomResidueName[1][i] << " not found " << endl;
754 0 : error ("missing Delta G parameter\n");
755 : }
756 : }
757 : }
758 0 : }
759 :
760 :
761 : //calculates neighbor list
762 24 : void SASA_HASEL::calcNlist() {
763 24 : if(!nopbc) makeWhole();
764 :
765 2904 : for(unsigned i = 0; i < natoms; i++) {
766 2880 : Nlist[i].clear();
767 : }
768 :
769 2904 : for(unsigned i = 0; i < natoms; i++) {
770 2880 : if (SASAparam[i].size()>0) {
771 103680 : for (unsigned j = 0; j < i; j++) {
772 101952 : if (SASAparam[j].size()>0) {
773 61344 : const Vector Delta_ij_vec = delta( getPosition(i), getPosition(j) );
774 61344 : double Delta_ij_mod = Delta_ij_vec.modulo()*10;
775 61344 : double overlapD = SASAparam[i][0]+SASAparam[j][0];
776 61344 : if (Delta_ij_mod < overlapD) {
777 26748 : Nlist.at(i).push_back (j);
778 26748 : Nlist.at(j).push_back (i);
779 : }
780 : }
781 : }
782 : }
783 : }
784 24 : }
785 :
786 :
787 : //calculates SASA according to Hasel et al., Tetrahedron Computer Methodology Vol. 1, No. 2, pp. 103-116, 1988
788 24 : void SASA_HASEL::calculate() {
789 24 : if(!nopbc) makeWhole();
790 :
791 24 : if(getExchangeStep()) nl_update = 0;
792 24 : if (firstStepFlag ==0) {
793 2 : readPDB();
794 2 : readSASAparam();
795 : }
796 24 : if (nl_update == 0) {
797 24 : calcNlist();
798 : }
799 :
800 :
801 24 : auto* moldat = plumed.getActionSet().selectLatest<GenericMolInfo*>(this);
802 24 : if( ! moldat ) error("Unable to find MOLINFO in input");
803 : double Si, sasai, bij;
804 24 : double sasa = 0;
805 24 : vector<Vector> derivatives( natoms );
806 2904 : for(unsigned i = 0; i < natoms; i++) {
807 2880 : derivatives[i][0] = 0.;
808 2880 : derivatives[i][1] = 0.;
809 2880 : derivatives[i][2] = 0.;
810 : }
811 :
812 24 : Tensor virial;
813 24 : vector <double> ddij_di(3);
814 24 : vector <double> dbij_di(3);
815 24 : vector <double> dAijt_di(3);
816 :
817 24 : if( sasa_type==TOTAL ) {
818 1452 : for(unsigned i = 0; i < natoms; i++) {
819 1440 : if(SASAparam[i].size() > 0) {
820 864 : double ri = SASAparam[i][0];
821 864 : Si = 4*M_PI*ri*ri;
822 : sasai = 1.0;
823 :
824 1728 : vector <vector <double> > derTerm( Nlist[i].size(), vector <double>(3));
825 :
826 864 : dAijt_di[0] = 0;
827 864 : dAijt_di[1] = 0;
828 864 : dAijt_di[2] = 0;
829 864 : int NumRes_i = moldat->getResidueNumber(atoms[i]);
830 :
831 27612 : for (unsigned j = 0; j < Nlist[i].size(); j++) {
832 : double pij = 0.3516;
833 :
834 26748 : int NumRes_j = moldat->getResidueNumber(atoms[Nlist[i][j]]);
835 26748 : if (NumRes_i==NumRes_j) {
836 4320 : if (CONNECTparam[i][0].compare(AtomResidueName[0][Nlist[i][j]])==0 || CONNECTparam[i][1].compare(AtomResidueName[0][Nlist[i][j]])==0 || CONNECTparam[i][2].compare(AtomResidueName[0][Nlist[i][j]])==0 || CONNECTparam[i][3].compare(AtomResidueName[0][Nlist[i][j]])==0) {
837 : pij = 0.8875;
838 : }
839 : }
840 26748 : if ( abs(NumRes_i-NumRes_j) == 1 ) {
841 10380 : if ((AtomResidueName[0][i] == "N" && AtomResidueName[0][Nlist[i][j]]== "CA") || (AtomResidueName[0][Nlist[i][j]] == "N" && AtomResidueName[0][i]== "CA")) {
842 : pij = 0.8875;
843 : }
844 : }
845 :
846 26748 : const Vector d_ij_vec = delta( getPosition(i), getPosition(Nlist[i][j]) );
847 26748 : double d_ij = d_ij_vec.modulo()*10;
848 :
849 26748 : double rj = SASAparam[Nlist[i][j]][0];
850 26748 : bij = M_PI*ri*(ri+rj-d_ij)*(1+(rj-ri)/d_ij); //Angstrom2
851 :
852 26748 : sasai = sasai*(1-SASAparam[i][1]*pij*bij/Si); //nondimensional
853 :
854 26748 : ddij_di[0] = -10*(getPosition(Nlist[i][j])[0]-getPosition(i)[0])/d_ij; //nondimensional
855 26748 : ddij_di[1] = -10*(getPosition(Nlist[i][j])[1]-getPosition(i)[1])/d_ij;
856 26748 : ddij_di[2] = -10*(getPosition(Nlist[i][j])[2]-getPosition(i)[2])/d_ij;
857 :
858 26748 : dbij_di[0] = -M_PI*ri*ddij_di[0]*(1+(ri+rj)*(rj-ri)/(d_ij*d_ij)); //Angstrom
859 26748 : dbij_di[1] = -M_PI*ri*ddij_di[1]*(1+(ri+rj)*(rj-ri)/(d_ij*d_ij));
860 26748 : dbij_di[2] = -M_PI*ri*ddij_di[2]*(1+(ri+rj)*(rj-ri)/(d_ij*d_ij));
861 :
862 26748 : dAijt_di[0] += -1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[0]; //Angstrom-1
863 26748 : dAijt_di[1] += -1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[1];
864 26748 : dAijt_di[2] += -1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[2];
865 :
866 26748 : derTerm[j][0] = 1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[0]; //Angstrom-1
867 26748 : derTerm[j][1] = 1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[1];
868 26748 : derTerm[j][2] = 1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[2];
869 :
870 : }
871 :
872 864 : sasa += Si*sasai/100; //nm2
873 :
874 864 : derivatives[i][0] += Si*sasai/10*dAijt_di[0]; //nm
875 864 : derivatives[i][1] += Si*sasai/10*dAijt_di[1];
876 864 : derivatives[i][2] += Si*sasai/10*dAijt_di[2];
877 :
878 27612 : for (unsigned j = 0; j < Nlist[i].size(); j++) {
879 26748 : derivatives[Nlist[i][j]][0] += Si*sasai/10*derTerm[j][0]; //nm
880 26748 : derivatives[Nlist[i][j]][1] += Si*sasai/10*derTerm[j][1];
881 26748 : derivatives[Nlist[i][j]][2] += Si*sasai/10*derTerm[j][2];
882 : }
883 864 : }
884 : }
885 : }
886 :
887 :
888 24 : if( sasa_type==TRANSFER ) {
889 :
890 12 : if (firstStepFlag ==0) {
891 1 : readMaxSurf();
892 : }
893 :
894 12 : if (firstStepFlag ==0 && DeltaGValues != "absent") {
895 1 : readDeltaG();
896 : }
897 :
898 12 : if (DeltaGValues == "absent") {
899 0 : computeDeltaG();
900 : }
901 :
902 :
903 1452 : for(unsigned i = 0; i < natoms; i++) {
904 :
905 :
906 :
907 1440 : if(SASAparam[i].size() > 0) {
908 864 : double ri = SASAparam[i][0];
909 864 : Si = 4*M_PI*ri*ri;
910 : sasai = 1.0;
911 :
912 1728 : vector <vector <double> > derTerm( Nlist[i].size(), vector <double>(3));
913 :
914 864 : dAijt_di[0] = 0;
915 864 : dAijt_di[1] = 0;
916 864 : dAijt_di[2] = 0;
917 864 : int NumRes_i = moldat->getResidueNumber(atoms[i]);
918 :
919 27612 : for (unsigned j = 0; j < Nlist[i].size(); j++) {
920 : double pij = 0.3516;
921 :
922 26748 : int NumRes_j = moldat->getResidueNumber(atoms[Nlist[i][j]]);
923 26748 : if (NumRes_i==NumRes_j) {
924 4320 : if (CONNECTparam[i][0].compare(AtomResidueName[0][Nlist[i][j]])==0 || CONNECTparam[i][1].compare(AtomResidueName[0][Nlist[i][j]])==0 || CONNECTparam[i][2].compare(AtomResidueName[0][Nlist[i][j]])==0 || CONNECTparam[i][3].compare(AtomResidueName[0][Nlist[i][j]])==0) {
925 : pij = 0.8875;
926 : }
927 : }
928 26748 : if ( abs(NumRes_i-NumRes_j) == 1 ) {
929 10380 : if ((AtomResidueName[0][i] == "N" && AtomResidueName[0][Nlist[i][j]]== "CA") || (AtomResidueName[0][Nlist[i][j]] == "N" && AtomResidueName[0][i]== "CA")) {
930 : pij = 0.8875;
931 : }
932 : }
933 :
934 26748 : const Vector d_ij_vec = delta( getPosition(i), getPosition(Nlist[i][j]) );
935 26748 : double d_ij = d_ij_vec.modulo()*10;
936 :
937 26748 : double rj = SASAparam[Nlist[i][j]][0];
938 26748 : bij = M_PI*ri*(ri+rj-d_ij)*(1+(rj-ri)/d_ij); //Angstrom2
939 :
940 26748 : sasai = sasai*(1-SASAparam[i][1]*pij*bij/Si); //nondimensional
941 :
942 26748 : ddij_di[0] = -10*(getPosition(Nlist[i][j])[0]-getPosition(i)[0])/d_ij; //nondimensional
943 26748 : ddij_di[1] = -10*(getPosition(Nlist[i][j])[1]-getPosition(i)[1])/d_ij;
944 26748 : ddij_di[2] = -10*(getPosition(Nlist[i][j])[2]-getPosition(i)[2])/d_ij;
945 :
946 26748 : dbij_di[0] = -M_PI*ri*ddij_di[0]*(1+(ri+rj)*(rj-ri)/(d_ij*d_ij)); //Angstrom
947 26748 : dbij_di[1] = -M_PI*ri*ddij_di[1]*(1+(ri+rj)*(rj-ri)/(d_ij*d_ij));
948 26748 : dbij_di[2] = -M_PI*ri*ddij_di[2]*(1+(ri+rj)*(rj-ri)/(d_ij*d_ij));
949 :
950 26748 : dAijt_di[0] += -1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[0]; //Angstrom-1
951 26748 : dAijt_di[1] += -1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[1];
952 26748 : dAijt_di[2] += -1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[2];
953 :
954 26748 : derTerm[j][0] = 1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[0]; //Angstrom-1
955 26748 : derTerm[j][1] = 1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[1];
956 26748 : derTerm[j][2] = 1/(Si/(SASAparam[i][1]*pij)-bij)*dbij_di[2];
957 :
958 : }
959 :
960 2880 : if (AtomResidueName[0][i] == "N" || AtomResidueName[0][i] == "CA" || AtomResidueName[0][i] == "C" || AtomResidueName[0][i] == "O" || AtomResidueName[0][i] == "H") {
961 :
962 720 : sasa += Si*sasai/MaxSurf[i][0]*DeltaG[natoms][0]; //kJ/mol
963 :
964 :
965 720 : derivatives[i][0] += Si*sasai*dAijt_di[0]/MaxSurf[i][0]*DeltaG[natoms][0]*10; //kJ/mol/nm
966 720 : derivatives[i][1] += Si*sasai*dAijt_di[1]/MaxSurf[i][0]*DeltaG[natoms][0]*10;
967 720 : derivatives[i][2] += Si*sasai*dAijt_di[2]/MaxSurf[i][0]*DeltaG[natoms][0]*10;
968 : }
969 :
970 2880 : if (AtomResidueName[0][i] != "N" && AtomResidueName[0][i] != "CA" && AtomResidueName[0][i] != "C" && AtomResidueName[0][i] != "O" && AtomResidueName[0][i] != "H") {
971 144 : sasa += Si*sasai/MaxSurf[i][1]*DeltaG[i][0]; //kJ/mol
972 :
973 144 : derivatives[i][0] += Si*sasai*dAijt_di[0]/MaxSurf[i][1]*DeltaG[i][0]*10; //kJ/mol/nm
974 144 : derivatives[i][1] += Si*sasai*dAijt_di[1]/MaxSurf[i][1]*DeltaG[i][0]*10;
975 144 : derivatives[i][2] += Si*sasai*dAijt_di[2]/MaxSurf[i][1]*DeltaG[i][0]*10;
976 : }
977 :
978 :
979 27612 : for (unsigned j = 0; j < Nlist[i].size(); j++) {
980 86883 : if (AtomResidueName[0][i] == "N" || AtomResidueName[0][i] == "CA" || AtomResidueName[0][i] == "C" || AtomResidueName[0][i] == "O" || AtomResidueName[0][i] == "H") {
981 22438 : derivatives[Nlist[i][j]][0] += Si*sasai*10*derTerm[j][0]/MaxSurf[i][0]*DeltaG[natoms][0]; //kJ/mol/nm
982 22438 : derivatives[Nlist[i][j]][1] += Si*sasai*10*derTerm[j][1]/MaxSurf[i][0]*DeltaG[natoms][0];
983 22438 : derivatives[Nlist[i][j]][2] += Si*sasai*10*derTerm[j][2]/MaxSurf[i][0]*DeltaG[natoms][0];
984 : }
985 :
986 86883 : if (AtomResidueName[0][i] != "N" && AtomResidueName[0][i] != "CA" && AtomResidueName[0][i] != "C" && AtomResidueName[0][i] != "O" && AtomResidueName[0][i] != "H") {
987 4310 : derivatives[Nlist[i][j]][0] += Si*sasai*10*derTerm[j][0]/MaxSurf[i][1]*DeltaG[i][0]; //kJ/mol/nm
988 4310 : derivatives[Nlist[i][j]][1] += Si*sasai*10*derTerm[j][1]/MaxSurf[i][1]*DeltaG[i][0];
989 4310 : derivatives[Nlist[i][j]][2] += Si*sasai*10*derTerm[j][2]/MaxSurf[i][1]*DeltaG[i][0];
990 : }
991 : }
992 864 : }
993 : }
994 : }
995 :
996 :
997 2904 : for(unsigned i=0; i<natoms; i++) {
998 2880 : setAtomsDerivatives(i,derivatives[i]);
999 2880 : virial -= Tensor(getPosition(i),derivatives[i]);
1000 : }
1001 :
1002 24 : setBoxDerivatives(virial);
1003 24 : setValue(sasa);
1004 24 : firstStepFlag = 1;
1005 24 : ++nl_update;
1006 24 : if (nl_update == stride) {
1007 24 : nl_update = 0;
1008 : }
1009 : // setBoxDerivativesNoPbc();
1010 24 : }
1011 :
1012 : }//namespace sasa
1013 : }//namespace PLMD
|