a6d2ede4a8d8e08c512cec66635c42d3b8d05be4
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / resall.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "resall.h"
40
41 #include <cctype>
42 #include <cstdlib>
43 #include <cstring>
44
45 #include <algorithm>
46 #include <string>
47 #include <vector>
48
49 #include "gromacs/gmxpreprocess/fflibutil.h"
50 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
51 #include "gromacs/gmxpreprocess/grompp_impl.h"
52 #include "gromacs/gmxpreprocess/notset.h"
53 #include "gromacs/gmxpreprocess/pgutil.h"
54 #include "gromacs/topology/atoms.h"
55 #include "gromacs/topology/symtab.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
60 #include "gromacs/utility/strdb.h"
61
62 #include "hackblock.h"
63
64 PreprocessingAtomTypes read_atype(const char* ffdir, t_symtab* tab)
65 {
66     FILE*   in;
67     char    buf[STRLEN], name[STRLEN];
68     double  m;
69     t_atom* a;
70
71     std::vector<std::string> files = fflib_search_file_end(ffdir, ".atp", TRUE);
72     snew(a, 1);
73     PreprocessingAtomTypes at;
74
75     for (const auto& filename : files)
76     {
77         in = fflib_open(filename);
78         while (!feof(in))
79         {
80             /* Skip blank or comment-only lines */
81             do
82             {
83                 if (fgets2(buf, STRLEN, in) != nullptr)
84                 {
85                     strip_comment(buf);
86                     trim(buf);
87                 }
88             } while ((feof(in) == 0) && strlen(buf) == 0);
89
90             if (sscanf(buf, "%s%lf", name, &m) == 2)
91             {
92                 a->m = m;
93                 at.addType(tab, *a, name, InteractionOfType({}, {}), 0, 0);
94                 fflush(stderr);
95             }
96             else
97             {
98                 gmx_fatal(FARGS, "Invalid atomtype format: '%s'", buf);
99             }
100         }
101         gmx_ffclose(in);
102     }
103     sfree(a);
104     return at;
105 }
106
107 static void print_resatoms(FILE* out, const PreprocessingAtomTypes& atype, const PreprocessResidue& rtpDBEntry)
108 {
109     fprintf(out, "[ %s ]\n", rtpDBEntry.resname.c_str());
110     fprintf(out, " [ atoms ]\n");
111
112     for (int j = 0; (j < rtpDBEntry.natom()); j++)
113     {
114         int         tp   = rtpDBEntry.atom[j].type;
115         const char* tpnm = atype.atomNameFromAtomType(tp);
116         if (tpnm == nullptr)
117         {
118             gmx_fatal(FARGS, "Incorrect atomtype (%d)", tp);
119         }
120         fprintf(out, "%6s  %6s  %8.3f  %6d\n", *(rtpDBEntry.atomname[j]), tpnm,
121                 rtpDBEntry.atom[j].q, rtpDBEntry.cgnr[j]);
122     }
123 }
124
125 static bool read_atoms(FILE* in, char* line, PreprocessResidue* r0, t_symtab* tab, PreprocessingAtomTypes* atype)
126 {
127     int    cg;
128     char   buf[256], buf1[256];
129     double q;
130
131     /* Read Atoms */
132     r0->atom.clear();
133     r0->atomname.clear();
134     r0->cgnr.clear();
135     while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
136     {
137         if (sscanf(line, "%s%s%lf%d", buf, buf1, &q, &cg) != 4)
138         {
139             return FALSE;
140         }
141         r0->atomname.push_back(put_symtab(tab, buf));
142         r0->atom.emplace_back();
143         r0->atom.back().q = q;
144         r0->cgnr.push_back(cg);
145         int j = atype->atomTypeFromName(buf1);
146         if (j == NOTSET)
147         {
148             gmx_fatal(FARGS,
149                       "Atom type %s (residue %s) not found in atomtype "
150                       "database",
151                       buf1, r0->resname.c_str());
152         }
153         r0->atom.back().type = j;
154         r0->atom.back().m    = atype->atomMassFromAtomType(j);
155     }
156
157     return TRUE;
158 }
159
160 static bool read_bondeds(int bt, FILE* in, char* line, PreprocessResidue* rtpDBEntry)
161 {
162     char str[STRLEN];
163
164     while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
165     {
166         int n = 0;
167         int ni;
168         rtpDBEntry->rb[bt].b.emplace_back();
169         BondedInteraction* newBond = &rtpDBEntry->rb[bt].b.back();
170         for (int j = 0; j < btsNiatoms[bt]; j++)
171         {
172             if (sscanf(line + n, "%s%n", str, &ni) == 1)
173             {
174                 newBond->a[j] = str;
175             }
176             else
177             {
178                 return FALSE;
179             }
180             n += ni;
181         }
182         while (isspace(line[n]))
183         {
184             n++;
185         }
186         rtrim(line + n);
187         newBond->s = line + n;
188     }
189
190     return TRUE;
191 }
192
193 static void print_resbondeds(FILE* out, int bt, const PreprocessResidue& rtpDBEntry)
194 {
195     if (!rtpDBEntry.rb[bt].b.empty())
196     {
197         fprintf(out, " [ %s ]\n", btsNames[bt]);
198
199         for (const auto& b : rtpDBEntry.rb[bt].b)
200         {
201             for (int j = 0; j < btsNiatoms[bt]; j++)
202             {
203                 fprintf(out, "%6s ", b.a[j].c_str());
204             }
205             if (!b.s.empty())
206             {
207                 fprintf(out, "    %s", b.s.c_str());
208             }
209             fprintf(out, "\n");
210         }
211     }
212 }
213
214 static void check_rtp(gmx::ArrayRef<const PreprocessResidue> rtpDBEntry, const std::string& libfn)
215 {
216     /* check for double entries, assuming list is already sorted */
217     for (auto it = rtpDBEntry.begin() + 1; it != rtpDBEntry.end(); it++)
218     {
219         auto prev = it - 1;
220         if (gmx::equalCaseInsensitive(prev->resname, it->resname))
221         {
222             fprintf(stderr, "WARNING double entry %s in file %s\n", it->resname.c_str(), libfn.c_str());
223         }
224     }
225 }
226
227 static int get_bt(char* header)
228 {
229     int i;
230
231     for (i = 0; i < ebtsNR; i++)
232     {
233         if (gmx_strcasecmp(btsNames[i], header) == 0)
234         {
235             return i;
236         }
237     }
238     return NOTSET;
239 }
240
241 /* print all the ebtsNR type numbers */
242 static void print_resall_header(FILE* out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
243 {
244     fprintf(out, "[ bondedtypes ]\n");
245     fprintf(out,
246             "; bonds  angles  dihedrals  impropers all_dihedrals nr_exclusions  HH14  "
247             "remove_dih\n");
248     fprintf(out, " %5d  %6d  %9d  %9d  %14d  %14d %14d %14d\n\n", rtpDBEntry[0].rb[0].type,
249             rtpDBEntry[0].rb[1].type, rtpDBEntry[0].rb[2].type, rtpDBEntry[0].rb[3].type,
250             static_cast<int>(rtpDBEntry[0].bKeepAllGeneratedDihedrals), rtpDBEntry[0].nrexcl,
251             static_cast<int>(rtpDBEntry[0].bGenerateHH14Interactions),
252             static_cast<int>(rtpDBEntry[0].bRemoveDihedralIfWithImproper));
253 }
254
255 void print_resall(FILE* out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry, const PreprocessingAtomTypes& atype)
256 {
257     if (rtpDBEntry.empty())
258     {
259         return;
260     }
261
262     print_resall_header(out, rtpDBEntry);
263
264     for (const auto& r : rtpDBEntry)
265     {
266         if (r.natom() > 0)
267         {
268             print_resatoms(out, atype, r);
269             for (int bt = 0; bt < ebtsNR; bt++)
270             {
271                 print_resbondeds(out, bt, r);
272             }
273         }
274     }
275 }
276
277 void readResidueDatabase(const std::string&              rrdb,
278                          std::vector<PreprocessResidue>* rtpDBEntry,
279                          PreprocessingAtomTypes*         atype,
280                          t_symtab*                       tab,
281                          bool                            bAllowOverrideRTP)
282 {
283     FILE* in;
284     char  filebase[STRLEN], line[STRLEN], header[STRLEN];
285     int   bt, nparam;
286     int   dum1, dum2, dum3;
287     bool  bNextResidue, bError;
288
289     fflib_filename_base(rrdb.c_str(), filebase, STRLEN);
290
291     in = fflib_open(rrdb);
292
293     PreprocessResidue header_settings;
294
295     /* these bonded parameters will overwritten be when  *
296      * there is a [ bondedtypes ] entry in the .rtp file */
297     header_settings.rb[ebtsBONDS].type  = 1; /* normal bonds     */
298     header_settings.rb[ebtsANGLES].type = 1; /* normal angles    */
299     header_settings.rb[ebtsPDIHS].type  = 1; /* normal dihedrals */
300     header_settings.rb[ebtsIDIHS].type  = 2; /* normal impropers */
301     header_settings.rb[ebtsEXCLS].type  = 1; /* normal exclusions */
302     header_settings.rb[ebtsCMAP].type   = 1; /* normal cmap torsions */
303
304     header_settings.bKeepAllGeneratedDihedrals    = FALSE;
305     header_settings.nrexcl                        = 3;
306     header_settings.bGenerateHH14Interactions     = TRUE;
307     header_settings.bRemoveDihedralIfWithImproper = TRUE;
308
309     /* Column 5 & 6 aren't really bonded types, but we include
310      * them here to avoid introducing a new section:
311      * Column 5 : This controls the generation of dihedrals from the bonding.
312      *            All possible dihedrals are generated automatically. A value of
313      *            1 here means that all these are retained. A value of
314      *            0 here requires generated dihedrals be removed if
315      *              * there are any dihedrals on the same central atoms
316      *                specified in the residue topology, or
317      *              * there are other identical generated dihedrals
318      *                sharing the same central atoms, or
319      *              * there are other generated dihedrals sharing the
320      *                same central bond that have fewer hydrogen atoms
321      * Column 6: Number of bonded neighbors to exclude.
322      * Column 7: Generate 1,4 interactions between two hydrogen atoms
323      * Column 8: Remove proper dihedrals if centered on the same bond
324      *           as an improper dihedral
325      */
326     get_a_line(in, line, STRLEN);
327     if (!get_header(line, header))
328     {
329         gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
330     }
331     if (gmx::equalCaseInsensitive("bondedtypes", header, 5))
332     {
333         get_a_line(in, line, STRLEN);
334         if ((nparam = sscanf(line, "%d %d %d %d %d %d %d %d", &header_settings.rb[ebtsBONDS].type,
335                              &header_settings.rb[ebtsANGLES].type,
336                              &header_settings.rb[ebtsPDIHS].type, &header_settings.rb[ebtsIDIHS].type,
337                              &dum1, &header_settings.nrexcl, &dum2, &dum3))
338             < 4)
339         {
340             gmx_fatal(FARGS, "need 4 to 8 parameters in the header of .rtp file %s at line:\n%s\n",
341                       rrdb.c_str(), line);
342         }
343         header_settings.bKeepAllGeneratedDihedrals    = (dum1 != 0);
344         header_settings.bGenerateHH14Interactions     = (dum2 != 0);
345         header_settings.bRemoveDihedralIfWithImproper = (dum3 != 0);
346         get_a_line(in, line, STRLEN);
347         if (nparam < 5)
348         {
349             fprintf(stderr, "Using default: not generating all possible dihedrals\n");
350             header_settings.bKeepAllGeneratedDihedrals = FALSE;
351         }
352         if (nparam < 6)
353         {
354             fprintf(stderr, "Using default: excluding 3 bonded neighbors\n");
355             header_settings.nrexcl = 3;
356         }
357         if (nparam < 7)
358         {
359             fprintf(stderr, "Using default: generating 1,4 H--H interactions\n");
360             header_settings.bGenerateHH14Interactions = TRUE;
361         }
362         if (nparam < 8)
363         {
364             fprintf(stderr,
365                     "Using default: removing proper dihedrals found on the same bond as a proper "
366                     "dihedral\n");
367             header_settings.bRemoveDihedralIfWithImproper = TRUE;
368         }
369     }
370     else
371     {
372         fprintf(stderr,
373                 "Reading .rtp file without '[ bondedtypes ]' directive,\n"
374                 "Will proceed as if the entry was:\n");
375         print_resall_header(stderr, gmx::arrayRefFromArray(&header_settings, 1));
376     }
377     /* We don't know the current size of rrtp, but simply realloc immediately */
378     auto oldArrayEnd = rtpDBEntry->end();
379     while (!feof(in))
380     {
381         /* Initialise rtp entry structure */
382         rtpDBEntry->push_back(header_settings);
383         PreprocessResidue* res = &rtpDBEntry->back();
384         if (!get_header(line, header))
385         {
386             gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
387         }
388         res->resname  = header;
389         res->filebase = filebase;
390
391         get_a_line(in, line, STRLEN);
392         bError       = FALSE;
393         bNextResidue = FALSE;
394         do
395         {
396             if (!get_header(line, header))
397             {
398                 bError = TRUE;
399             }
400             else
401             {
402                 bt = get_bt(header);
403                 if (bt != NOTSET)
404                 {
405                     /* header is an bonded directive */
406                     bError = !read_bondeds(bt, in, line, res);
407                 }
408                 else if (gmx::equalCaseInsensitive("atoms", header, 5))
409                 {
410                     /* header is the atoms directive */
411                     bError = !read_atoms(in, line, res, tab, atype);
412                 }
413                 else
414                 {
415                     /* else header must be a residue name */
416                     bNextResidue = TRUE;
417                 }
418             }
419             if (bError)
420             {
421                 gmx_fatal(FARGS, "in .rtp file in residue %s at line:\n%s\n", res->resname.c_str(), line);
422             }
423         } while ((feof(in) == 0) && !bNextResidue);
424
425         if (res->natom() == 0)
426         {
427             gmx_fatal(FARGS, "No atoms found in .rtp file in residue %s\n", res->resname.c_str());
428         }
429
430         auto found = std::find_if(rtpDBEntry->begin(), rtpDBEntry->end() - 1,
431                                   [&res](const PreprocessResidue& entry) {
432                                       return gmx::equalCaseInsensitive(entry.resname, res->resname);
433                                   });
434
435         if (found != rtpDBEntry->end() - 1)
436         {
437             if (found >= oldArrayEnd)
438             {
439                 gmx_fatal(FARGS, "Found a second entry for '%s' in '%s'", res->resname.c_str(),
440                           rrdb.c_str());
441             }
442             if (bAllowOverrideRTP)
443             {
444                 fprintf(stderr,
445                         "Found another rtp entry for '%s' in '%s', ignoring this entry and keeping "
446                         "the one from '%s.rtp'\n",
447                         res->resname.c_str(), rrdb.c_str(), found->filebase.c_str());
448                 /* We should free all the data for this entry.
449                  * The current code gives a lot of dangling pointers.
450                  */
451                 rtpDBEntry->erase(rtpDBEntry->end() - 1);
452             }
453             else
454             {
455                 gmx_fatal(FARGS,
456                           "Found rtp entries for '%s' in both '%s' and '%s'. If you want the first "
457                           "definition to override the second one, set the -rtpo option of pdb2gmx.",
458                           res->resname.c_str(), found->filebase.c_str(), rrdb.c_str());
459             }
460         }
461     }
462     gmx_ffclose(in);
463
464     std::sort(rtpDBEntry->begin(), rtpDBEntry->end(), [](const PreprocessResidue& a, const PreprocessResidue& b) {
465         return std::lexicographical_compare(
466                 a.resname.begin(), a.resname.end(), b.resname.begin(), b.resname.end(),
467                 [](const char& c1, const char& c2) { return std::toupper(c1) < std::toupper(c2); });
468     });
469
470     check_rtp(*rtpDBEntry, rrdb);
471 }
472
473 /************************************************************
474  *
475  *                  SEARCH   ROUTINES
476  *
477  ***********************************************************/
478 static bool is_sign(char c)
479 {
480     return (c == '+' || c == '-');
481 }
482
483 /* Compares if the strins match except for a sign at the end
484  * of (only) one of the two.
485  */
486 static int neq_str_sign(const char* a1, const char* a2)
487 {
488     int l1, l2, lm;
489
490     l1 = static_cast<int>(strlen(a1));
491     l2 = static_cast<int>(strlen(a2));
492     lm = std::min(l1, l2);
493
494     if (lm >= 1 && ((l1 == l2 + 1 && is_sign(a1[l1 - 1])) || (l2 == l1 + 1 && is_sign(a2[l2 - 1])))
495         && gmx::equalCaseInsensitive(a1, a2, lm))
496     {
497         return lm;
498     }
499     else
500     {
501         return 0;
502     }
503 }
504
505 std::string searchResidueDatabase(const std::string& key, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
506 {
507     int         nbest, best, besti;
508     std::string bestbuf;
509
510     nbest = 0;
511     besti = -1;
512     /* We want to match at least one character */
513     best = 1;
514     for (auto it = rtpDBEntry.begin(); it != rtpDBEntry.end(); it++)
515     {
516         if (gmx::equalCaseInsensitive(key, it->resname))
517         {
518             besti = std::distance(rtpDBEntry.begin(), it);
519             nbest = 1;
520             break;
521         }
522         else
523         {
524             /* Allow a mismatch of at most a sign character (with warning) */
525             int n = neq_str_sign(key.c_str(), it->resname.c_str());
526             if (n >= best && n + 1 >= gmx::index(key.length()) && n + 1 >= gmx::index(it->resname.length()))
527             {
528                 if (n == best)
529                 {
530                     if (nbest == 1)
531                     {
532                         bestbuf = rtpDBEntry[besti].resname;
533                     }
534                     if (nbest >= 1)
535                     {
536                         bestbuf.append(" or ");
537                         bestbuf.append(it->resname);
538                     }
539                 }
540                 else
541                 {
542                     nbest = 0;
543                 }
544                 besti = std::distance(rtpDBEntry.begin(), it);
545                 best  = n;
546                 nbest++;
547             }
548         }
549     }
550     if (nbest > 1)
551     {
552         gmx_fatal(FARGS, "Residue '%s' not found in residue topology database, looks a bit like %s",
553                   key.c_str(), bestbuf.c_str());
554     }
555     else if (besti == -1)
556     {
557         gmx_fatal(FARGS, "Residue '%s' not found in residue topology database", key.c_str());
558     }
559     if (!gmx::equalCaseInsensitive(rtpDBEntry[besti].resname, key))
560     {
561         fprintf(stderr,
562                 "\nWARNING: '%s' not found in residue topology database, "
563                 "trying to use '%s'\n\n",
564                 key.c_str(), rtpDBEntry[besti].resname.c_str());
565     }
566
567     return rtpDBEntry[besti].resname;
568 }
569
570 gmx::ArrayRef<const PreprocessResidue>::const_iterator
571 getDatabaseEntry(const std::string& rtpname, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
572 {
573     auto found = std::find_if(rtpDBEntry.begin(), rtpDBEntry.end(),
574                               [&rtpname](const PreprocessResidue& entry) {
575                                   return gmx::equalCaseInsensitive(rtpname, entry.resname);
576                               });
577     if (found == rtpDBEntry.end())
578     {
579         /* This should never happen, since searchResidueDatabase should have been called
580          * before calling getDatabaseEntry.
581          */
582         gmx_fatal(FARGS, "Residue type '%s' not found in residue topology database", rtpname.c_str());
583     }
584
585     return found;
586 }