fac51b47a22947911ab3df0c0b24035a819a98c9
[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     int                      nratt = 0;
70     t_atom                  *a;
71
72     std::vector<std::string> files = fflib_search_file_end(ffdir, ".atp", TRUE);
73     snew(a, 1);
74     PreprocessingAtomTypes   at;
75
76     for (const auto &filename : files)
77     {
78         in = fflib_open(filename);
79         while (!feof(in))
80         {
81             /* Skip blank or comment-only lines */
82             do
83             {
84                 if (fgets2(buf, STRLEN, in) != nullptr)
85                 {
86                     strip_comment(buf);
87                     trim(buf);
88                 }
89             }
90             while ((feof(in) == 0) && strlen(buf) == 0);
91
92             if (sscanf(buf, "%s%lf", name, &m) == 2)
93             {
94                 a->m = m;
95                 at.addType(tab, *a, name, InteractionType({}, {}), 0, 0);
96                 fprintf(stderr, "\rAtomtype %d", ++nratt);
97                 fflush(stderr);
98             }
99             else
100             {
101                 fprintf(stderr, "\nInvalid format: %s\n", buf);
102             }
103         }
104         gmx_ffclose(in);
105     }
106     fprintf(stderr, "\n");
107     sfree(a);
108     return at;
109 }
110
111 static void print_resatoms(FILE *out, const PreprocessingAtomTypes &atype, const PreprocessResidue &rtpDBEntry)
112 {
113     fprintf(out, "[ %s ]\n", rtpDBEntry.resname.c_str());
114     fprintf(out, " [ atoms ]\n");
115
116     for (int j = 0; (j < rtpDBEntry.natom()); j++)
117     {
118         int         tp   = rtpDBEntry.atom[j].type;
119         const char *tpnm = atype.atomNameFromAtomType(tp);
120         if (tpnm == nullptr)
121         {
122             gmx_fatal(FARGS, "Incorrect atomtype (%d)", tp);
123         }
124         fprintf(out, "%6s  %6s  %8.3f  %6d\n",
125                 *(rtpDBEntry.atomname[j]), tpnm, rtpDBEntry.atom[j].q, rtpDBEntry.cgnr[j]);
126     }
127 }
128
129 static bool read_atoms(FILE *in, char *line,
130                        PreprocessResidue *r0, t_symtab *tab, PreprocessingAtomTypes *atype)
131 {
132     int    cg;
133     char   buf[256], buf1[256];
134     double q;
135
136     /* Read Atoms */
137     r0->atom.clear();
138     r0->atomname.clear();
139     r0->cgnr.clear();
140     while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
141     {
142         if (sscanf(line, "%s%s%lf%d", buf, buf1, &q, &cg) != 4)
143         {
144             return FALSE;
145         }
146         r0->atomname.push_back(put_symtab(tab, buf));
147         r0->atom.emplace_back();
148         r0->atom.back().q = q;
149         r0->cgnr.push_back(cg);
150         int j               = atype->atomTypeFromName(buf1);
151         if (j == NOTSET)
152         {
153             gmx_fatal(FARGS, "Atom type %s (residue %s) not found in atomtype "
154                       "database", buf1, r0->resname.c_str());
155         }
156         r0->atom.back().type = j;
157         r0->atom.back().m    = atype->atomMassAFromAtomType(j);
158     }
159
160     return TRUE;
161 }
162
163 static bool read_bondeds(int bt, FILE *in, char *line, PreprocessResidue *rtpDBEntry)
164 {
165     char str[STRLEN];
166
167     while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
168     {
169         int                n = 0;
170         int                ni;
171         rtpDBEntry->rb[bt].b.emplace_back();
172         BondedInteraction *newBond = &rtpDBEntry->rb[bt].b.back();
173         for (int j = 0; j < btsNiatoms[bt]; j++)
174         {
175             if (sscanf(line+n, "%s%n", str, &ni) == 1)
176             {
177                 newBond->a[j]  = str;
178             }
179             else
180             {
181                 return FALSE;
182             }
183             n += ni;
184         }
185         while (isspace(line[n]))
186         {
187             n++;
188         }
189         rtrim(line+n);
190         newBond->s = line+n;
191     }
192
193     return TRUE;
194 }
195
196 static void print_resbondeds(FILE *out, int bt, const PreprocessResidue &rtpDBEntry)
197 {
198     if (!rtpDBEntry.rb[bt].b.empty())
199     {
200         fprintf(out, " [ %s ]\n", btsNames[bt]);
201
202         for (const auto &b : rtpDBEntry.rb[bt].b)
203         {
204             for (int j = 0; j < btsNiatoms[bt]; j++)
205             {
206                 fprintf(out, "%6s ", b.a[j].c_str());
207             }
208             if (!b.s.empty())
209             {
210                 fprintf(out, "    %s", b.s.c_str());
211             }
212             fprintf(out, "\n");
213         }
214     }
215 }
216
217 static void check_rtp(gmx::ArrayRef<const PreprocessResidue> rtpDBEntry, const std::string &libfn)
218 {
219     /* check for double entries, assuming list is already sorted */
220     for (auto it = rtpDBEntry.begin() + 1; it != rtpDBEntry.end(); it++)
221     {
222         auto prev = it-1;
223         if (gmx::equalCaseInsensitive(prev->resname, it->resname))
224         {
225             fprintf(stderr, "WARNING double entry %s in file %s\n",
226                     it->resname.c_str(), libfn.c_str());
227         }
228     }
229 }
230
231 static int get_bt(char* header)
232 {
233     int i;
234
235     for (i = 0; i < ebtsNR; i++)
236     {
237         if (gmx_strcasecmp(btsNames[i], header) == 0)
238         {
239             return i;
240         }
241     }
242     return NOTSET;
243 }
244
245 /* print all the ebtsNR type numbers */
246 static void print_resall_header(FILE *out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
247 {
248     fprintf(out, "[ bondedtypes ]\n");
249     fprintf(out, "; bonds  angles  dihedrals  impropers all_dihedrals nr_exclusions  HH14  remove_dih\n");
250     fprintf(out, " %5d  %6d  %9d  %9d  %14d  %14d %14d %14d\n\n",
251             rtpDBEntry[0].rb[0].type,
252             rtpDBEntry[0].rb[1].type,
253             rtpDBEntry[0].rb[2].type,
254             rtpDBEntry[0].rb[3].type,
255             static_cast<int>(rtpDBEntry[0].bKeepAllGeneratedDihedrals),
256             rtpDBEntry[0].nrexcl,
257             static_cast<int>(rtpDBEntry[0].bGenerateHH14Interactions),
258             static_cast<int>(rtpDBEntry[0].bRemoveDihedralIfWithImproper));
259 }
260
261 void print_resall(FILE *out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry,
262                   const PreprocessingAtomTypes &atype)
263 {
264     if (rtpDBEntry.empty())
265     {
266         return;
267     }
268
269     print_resall_header(out, rtpDBEntry);
270
271     for (const auto &r : rtpDBEntry)
272     {
273         if (r.natom() > 0)
274         {
275             print_resatoms(out, atype, r);
276             for (int bt = 0; bt < ebtsNR; bt++)
277             {
278                 print_resbondeds(out, bt, r);
279             }
280         }
281     }
282 }
283
284 void readResidueDatabase(const std::string &rrdb, std::vector<PreprocessResidue> *rtpDBEntry,
285                          PreprocessingAtomTypes *atype, t_symtab *tab,
286                          bool bAllowOverrideRTP)
287 {
288     FILE         *in;
289     char          filebase[STRLEN], line[STRLEN], header[STRLEN];
290     int           bt, nparam;
291     int           dum1, dum2, dum3;
292     bool          bNextResidue, bError;
293
294     fflib_filename_base(rrdb.c_str(), filebase, STRLEN);
295
296     in = fflib_open(rrdb);
297
298     PreprocessResidue header_settings;
299
300     /* these bonded parameters will overwritten be when  *
301      * there is a [ bondedtypes ] entry in the .rtp file */
302     header_settings.rb[ebtsBONDS].type  = 1; /* normal bonds     */
303     header_settings.rb[ebtsANGLES].type = 1; /* normal angles    */
304     header_settings.rb[ebtsPDIHS].type  = 1; /* normal dihedrals */
305     header_settings.rb[ebtsIDIHS].type  = 2; /* normal impropers */
306     header_settings.rb[ebtsEXCLS].type  = 1; /* normal exclusions */
307     header_settings.rb[ebtsCMAP].type   = 1; /* normal cmap torsions */
308
309     header_settings.bKeepAllGeneratedDihedrals    = FALSE;
310     header_settings.nrexcl                        = 3;
311     header_settings.bGenerateHH14Interactions     = TRUE;
312     header_settings.bRemoveDihedralIfWithImproper = TRUE;
313
314     /* Column 5 & 6 aren't really bonded types, but we include
315      * them here to avoid introducing a new section:
316      * Column 5 : This controls the generation of dihedrals from the bonding.
317      *            All possible dihedrals are generated automatically. A value of
318      *            1 here means that all these are retained. A value of
319      *            0 here requires generated dihedrals be removed if
320      *              * there are any dihedrals on the same central atoms
321      *                specified in the residue topology, or
322      *              * there are other identical generated dihedrals
323      *                sharing the same central atoms, or
324      *              * there are other generated dihedrals sharing the
325      *                same central bond that have fewer hydrogen atoms
326      * Column 6: Number of bonded neighbors to exclude.
327      * Column 7: Generate 1,4 interactions between two hydrogen atoms
328      * Column 8: Remove proper dihedrals if centered on the same bond
329      *           as an improper dihedral
330      */
331     get_a_line(in, line, STRLEN);
332     if (!get_header(line, header))
333     {
334         gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
335     }
336     if (gmx_strncasecmp("bondedtypes", header, 5) == 0)
337     {
338         get_a_line(in, line, STRLEN);
339         if ((nparam = sscanf(line, "%d %d %d %d %d %d %d %d",
340                              &header_settings.rb[ebtsBONDS].type, &header_settings.rb[ebtsANGLES].type,
341                              &header_settings.rb[ebtsPDIHS].type, &header_settings.rb[ebtsIDIHS].type,
342                              &dum1, &header_settings.nrexcl, &dum2, &dum3)) < 4)
343         {
344             gmx_fatal(FARGS, "need 4 to 8 parameters in the header of .rtp file %s at line:\n%s\n", rrdb.c_str(), line);
345         }
346         header_settings.bKeepAllGeneratedDihedrals    = (dum1 != 0);
347         header_settings.bGenerateHH14Interactions     = (dum2 != 0);
348         header_settings.bRemoveDihedralIfWithImproper = (dum3 != 0);
349         get_a_line(in, line, STRLEN);
350         if (nparam < 5)
351         {
352             fprintf(stderr, "Using default: not generating all possible dihedrals\n");
353             header_settings.bKeepAllGeneratedDihedrals = FALSE;
354         }
355         if (nparam < 6)
356         {
357             fprintf(stderr, "Using default: excluding 3 bonded neighbors\n");
358             header_settings.nrexcl = 3;
359         }
360         if (nparam < 7)
361         {
362             fprintf(stderr, "Using default: generating 1,4 H--H interactions\n");
363             header_settings.bGenerateHH14Interactions = TRUE;
364         }
365         if (nparam < 8)
366         {
367             fprintf(stderr, "Using default: removing proper dihedrals found on the same bond as a proper dihedral\n");
368             header_settings.bRemoveDihedralIfWithImproper = TRUE;
369         }
370     }
371     else
372     {
373         fprintf(stderr,
374                 "Reading .rtp file without '[ bondedtypes ]' directive,\n"
375                 "Will proceed as if the entry was:\n");
376         print_resall_header(stderr, gmx::arrayRefFromArray(&header_settings, 1));
377     }
378     /* We don't know the current size of rrtp, but simply realloc immediately */
379     auto oldArrayEnd = rtpDBEntry->end() - 1;
380     while (!feof(in))
381     {
382         /* Initialise rtp entry structure */
383         rtpDBEntry->push_back(header_settings);
384         PreprocessResidue *res = &rtpDBEntry->back();
385         if (!get_header(line, header))
386         {
387             gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
388         }
389         res->resname  = header;
390         res->filebase = filebase;
391
392         get_a_line(in, line, STRLEN);
393         bError       = FALSE;
394         bNextResidue = FALSE;
395         do
396         {
397             if (!get_header(line, header))
398             {
399                 bError = TRUE;
400             }
401             else
402             {
403                 bt = get_bt(header);
404                 if (bt != NOTSET)
405                 {
406                     /* header is an bonded directive */
407                     bError = !read_bondeds(bt, in, line, res);
408                 }
409                 else if (gmx_strncasecmp("atoms", header, 5) == 0)
410                 {
411                     /* header is the atoms directive */
412                     bError = !read_atoms(in, line, res, tab, atype);
413                 }
414                 else
415                 {
416                     /* else header must be a residue name */
417                     bNextResidue = TRUE;
418                 }
419             }
420             if (bError)
421             {
422                 gmx_fatal(FARGS, "in .rtp file in residue %s at line:\n%s\n",
423                           res->resname.c_str(), line);
424             }
425         }
426         while ((feof(in) == 0) && !bNextResidue);
427
428         if (res->natom() == 0)
429         {
430             gmx_fatal(FARGS, "No atoms found in .rtp file in residue %s\n",
431                       res->resname.c_str());
432         }
433
434         auto found = std::find_if(rtpDBEntry->begin(), rtpDBEntry->end()-1,
435                                   [&res](const PreprocessResidue &entry)
436                                   { return gmx::equalCaseInsensitive(entry.resname, res->resname); });
437
438         if (found == rtpDBEntry->end() - 1)
439         {
440             int size = rtpDBEntry->size() - 1;
441             fprintf(stderr, "\rResidue %d", size);
442             fflush(stderr);
443         }
444         else
445         {
446             if (found >= oldArrayEnd)
447             {
448                 gmx_fatal(FARGS, "Found a second entry for '%s' in '%s'",
449                           res->resname.c_str(), rrdb.c_str());
450             }
451             if (bAllowOverrideRTP)
452             {
453                 fprintf(stderr, "\n");
454                 fprintf(stderr, "Found another rtp entry for '%s' in '%s', ignoring this entry and keeping the one from '%s.rtp'\n",
455                         res->resname.c_str(), rrdb.c_str(), found->filebase.c_str());
456                 /* We should free all the data for this entry.
457                  * The current code gives a lot of dangling pointers.
458                  */
459                 rtpDBEntry->erase(rtpDBEntry->end() - 1);
460             }
461             else
462             {
463                 gmx_fatal(FARGS, "Found rtp entries for '%s' in both '%s' and '%s'. If you want the first definition to override the second one, set the -rtpo option of pdb2gmx.", res->resname.c_str(), found->filebase.c_str(), rrdb.c_str());
464             }
465         }
466     }
467     gmx_ffclose(in);
468
469     fprintf(stderr, "\nSorting it all out...\n");
470     std::sort(rtpDBEntry->begin(), rtpDBEntry->end(), []
471                   (const PreprocessResidue &a,
472                   const PreprocessResidue &b)
473               {return std::lexicographical_compare(
474                        a.resname.begin(), a.resname.end(),
475                        b.resname.begin(), b.resname.end(),
476                        [](const char &c1, const char &c2)
477                        { return std::toupper(c1) < std::toupper(c2); }); });
478
479     check_rtp(*rtpDBEntry, rrdb);
480 }
481
482 /************************************************************
483  *
484  *                  SEARCH   ROUTINES
485  *
486  ***********************************************************/
487 static bool is_sign(char c)
488 {
489     return (c == '+' || c == '-');
490 }
491
492 /* Compares if the strins match except for a sign at the end
493  * of (only) one of the two.
494  */
495 static int neq_str_sign(const char *a1, const char *a2)
496 {
497     int l1, l2, lm;
498
499     l1 = static_cast<int>(strlen(a1));
500     l2 = static_cast<int>(strlen(a2));
501     lm = std::min(l1, l2);
502
503     if (lm >= 1 &&
504         ((l1 == l2+1 && is_sign(a1[l1-1])) ||
505          (l2 == l1+1 && is_sign(a2[l2-1]))) &&
506         gmx_strncasecmp(a1, a2, lm) == 0)
507     {
508         return lm;
509     }
510     else
511     {
512         return 0;
513     }
514 }
515
516 std::string searchResidueDatabase(const std::string &key, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
517 {
518     int         nbest, best, besti;
519     std::string bestbuf;
520
521     nbest =  0;
522     besti = -1;
523     /* We want to match at least one character */
524     best  =  1;
525     for (auto it = rtpDBEntry.begin(); it != rtpDBEntry.end(); it++)
526     {
527         if (gmx::equalCaseInsensitive(key, it->resname))
528         {
529             besti = std::distance(rtpDBEntry.begin(), it);
530             nbest = 1;
531             break;
532         }
533         else
534         {
535             /* Allow a mismatch of at most a sign character (with warning) */
536             int n = neq_str_sign(key.c_str(), it->resname.c_str());
537             if (n >= best &&
538                 n+1 >= gmx::index(key.length()) &&
539                 n+1 >= gmx::index(it->resname.length()))
540             {
541                 if (n == best)
542                 {
543                     if (nbest == 1)
544                     {
545                         bestbuf = rtpDBEntry[besti].resname;
546                     }
547                     if (nbest >= 1)
548                     {
549                         bestbuf.append(" or ");
550                         bestbuf.append(it->resname);
551                     }
552                 }
553                 else
554                 {
555                     nbest = 0;
556                 }
557                 besti = std::distance(rtpDBEntry.begin(), it);
558                 best  = n;
559                 nbest++;
560             }
561         }
562     }
563     if (nbest > 1)
564     {
565         gmx_fatal(FARGS, "Residue '%s' not found in residue topology database, looks a bit like %s", key.c_str(), bestbuf.c_str());
566     }
567     else if (besti == -1)
568     {
569         gmx_fatal(FARGS, "Residue '%s' not found in residue topology database", key.c_str());
570     }
571     if (!gmx::equalCaseInsensitive(rtpDBEntry[besti].resname, key))
572     {
573         fprintf(stderr,
574                 "\nWARNING: '%s' not found in residue topology database, "
575                 "trying to use '%s'\n\n", key.c_str(), rtpDBEntry[besti].resname.c_str());
576     }
577
578     return rtpDBEntry[besti].resname;
579 }
580
581 gmx::ArrayRef<const PreprocessResidue>::const_iterator
582 getDatabaseEntry(const std::string &rtpname, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
583 {
584     auto found = std::find_if(rtpDBEntry.begin(), rtpDBEntry.end(),
585                               [&rtpname](const PreprocessResidue &entry)
586                               { return gmx::equalCaseInsensitive(rtpname, entry.resname); });
587     if (found == rtpDBEntry.end())
588     {
589         /* This should never happen, since searchResidueDatabase should have been called
590          * before calling getDatabaseEntry.
591          */
592         gmx_fatal(FARGS, "Residue type '%s' not found in residue topology database", rtpname.c_str());
593     }
594
595     return found;
596 }