Remove unused thole polarization rfac parameter
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / ter_db.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,2017,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "ter_db.h"
41
42 #include <array>
43 #include <cctype>
44 #include <cstring>
45
46 #include <algorithm>
47 #include <optional>
48 #include <string>
49 #include <vector>
50
51 #include "gromacs/fileio/gmxfio.h"
52 #include "gromacs/gmxpreprocess/fflibutil.h"
53 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
54 #include "gromacs/gmxpreprocess/grompp_impl.h"
55 #include "gromacs/gmxpreprocess/h_db.h"
56 #include "gromacs/gmxpreprocess/notset.h"
57 #include "gromacs/gmxpreprocess/toputil.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/enumerationhelpers.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
63 #include "gromacs/utility/strdb.h"
64 #include "gromacs/utility/stringtoenumvalueconverter.h"
65
66 #include "hackblock.h"
67 #include "resall.h"
68
69 /* use bonded types definitions in hackblock.h */
70 enum class ReplaceType : int
71 {
72     Repl,
73     Add,
74     Del,
75     Count
76 };
77
78 static const char* enumValueToString(ReplaceType enumValue)
79 {
80     constexpr gmx::EnumerationArray<ReplaceType, const char*> replaceTypeNames = { "replace",
81                                                                                    "add",
82                                                                                    "delete" };
83     return replaceTypeNames[enumValue];
84 }
85
86 template<typename EnumType>
87 static std::optional<EnumType> findTypeFromKeyword(char* keyw)
88 {
89     gmx::StringToEnumValueConverter<EnumType, enumValueToString, gmx::StringCompareType::CaseInsensitive, gmx::StripStrings::Yes> converter;
90     return converter.valueFrom(keyw);
91 }
92
93 static void read_atom(char* line, bool bAdd, std::string* nname, t_atom* a, PreprocessingAtomTypes* atype, int* cgnr)
94 {
95     int    nr, i;
96     char   buf[5][30];
97     double m, q;
98
99     /* This code is messy, because of support for different formats:
100      * for replace: [new name] <atom type> <m> <q> [cgnr (old format)]
101      * for add:                <atom type> <m> <q> [cgnr]
102      */
103     nr = sscanf(line, "%s %s %s %s %s", buf[0], buf[1], buf[2], buf[3], buf[4]);
104
105     /* Here there an ambiguity due to the old replace format with cgnr,
106      * which was read for years, but ignored in the rest of the code.
107      * We have to assume that the atom type does not start with a digit
108      * to make a line with 4 entries uniquely interpretable.
109      */
110     if (!bAdd && nr == 4 && isdigit(buf[1][0]))
111     {
112         nr = 3;
113     }
114
115     if (nr < 3 || nr > 4)
116     {
117         gmx_fatal(FARGS,
118                   "Reading Termini Database: expected %d or %d items of atom data in stead of %d "
119                   "on line\n%s",
120                   3,
121                   4,
122                   nr,
123                   line);
124     }
125     i = 0;
126     if (!bAdd)
127     {
128         if (nr == 4)
129         {
130             *nname = buf[i++];
131         }
132         else
133         {
134             *nname = "";
135         }
136     }
137     a->type = *atype->atomTypeFromName(buf[i++]);
138     sscanf(buf[i++], "%lf", &m);
139     a->m = m;
140     sscanf(buf[i++], "%lf", &q);
141     a->q = q;
142     if (bAdd && nr == 4)
143     {
144         sscanf(buf[i++], "%d", cgnr);
145     }
146     else
147     {
148         *cgnr = NOTSET;
149     }
150 }
151
152 static void print_atom(FILE* out, const t_atom& a, PreprocessingAtomTypes* atype)
153 {
154     fprintf(out, "\t%s\t%g\t%g\n", *atype->atomNameFromAtomType(a.type), a.m, a.q);
155 }
156
157 static void print_ter_db(const char*                                ff,
158                          char                                       C,
159                          gmx::ArrayRef<const MoleculePatchDatabase> tb,
160                          PreprocessingAtomTypes*                    atype)
161 {
162     std::string buf = gmx::formatString("%s-%c.tdb", ff, C);
163     FILE*       out = gmx_fio_fopen(buf.c_str(), "w");
164
165     for (const auto& modification : tb)
166     {
167         fprintf(out, "[ %s ]\n", modification.name.c_str());
168
169         if (std::any_of(modification.hack.begin(), modification.hack.end(), [](const auto& mod) {
170                 return mod.type() == MoleculePatchType::Replace;
171             }))
172         {
173             fprintf(out, "[ %s ]\n", enumValueToString(ReplaceType::Repl));
174             for (const auto& hack : modification.hack)
175             {
176                 if (hack.type() == MoleculePatchType::Replace)
177                 {
178                     fprintf(out, "%s\t", hack.oname.c_str());
179                     print_atom(out, hack.atom.back(), atype);
180                 }
181             }
182         }
183         if (std::any_of(modification.hack.begin(), modification.hack.end(), [](const auto& mod) {
184                 return mod.type() == MoleculePatchType::Add;
185             }))
186         {
187             fprintf(out, "[ %s ]\n", enumValueToString(ReplaceType::Add));
188             for (const auto& hack : modification.hack)
189             {
190                 if (hack.type() == MoleculePatchType::Add)
191                 {
192                     print_ab(out, hack, hack.nname.c_str());
193                     print_atom(out, hack.atom.back(), atype);
194                 }
195             }
196         }
197         if (std::any_of(modification.hack.begin(), modification.hack.end(), [](const auto& mod) {
198                 return mod.type() == MoleculePatchType::Delete;
199             }))
200         {
201             fprintf(out, "[ %s ]\n", enumValueToString(ReplaceType::Del));
202             for (const auto& hack : modification.hack)
203             {
204                 if (hack.type() == MoleculePatchType::Delete)
205                 {
206                     fprintf(out, "%s\n", hack.oname.c_str());
207                 }
208             }
209         }
210         for (auto bt : gmx::EnumerationWrapper<BondedTypes>{})
211         {
212             if (!modification.rb[bt].b.empty())
213             {
214                 fprintf(out, "[ %s ]\n", enumValueToString(bt));
215                 for (const auto& b : modification.rb[bt].b)
216                 {
217                     for (int k = 0; k < enumValueToNumIAtoms(bt); k++)
218                     {
219                         fprintf(out, "%s%s", k ? "\t" : "", b.a[k].c_str());
220                     }
221                     if (!b.s.empty())
222                     {
223                         fprintf(out, "\t%s", b.s.c_str());
224                     }
225                     fprintf(out, "\n");
226                 }
227             }
228         }
229         fprintf(out, "\n");
230     }
231     gmx_fio_fclose(out);
232 }
233
234 static void read_ter_db_file(const char*                         fn,
235                              std::vector<MoleculePatchDatabase>* tbptr,
236                              PreprocessingAtomTypes*             atype)
237 {
238     char filebase[STRLEN];
239     char header[STRLEN], buf[STRLEN], line[STRLEN];
240
241     fflib_filename_base(fn, filebase, STRLEN);
242     /* Remove the C/N termini extension */
243     char* ptr = strrchr(filebase, '.');
244     if (ptr != nullptr)
245     {
246         ptr[0] = '\0';
247     }
248
249     FILE* in = fflib_open(fn);
250
251     std::optional<BondedTypes> btkw;
252     std::optional<ReplaceType> rtkw;
253     get_a_line(in, line, STRLEN);
254     MoleculePatchDatabase* block = nullptr;
255     while (!feof(in))
256     {
257         if (get_header(line, header))
258         {
259             /* this is a new block, or a new keyword */
260             btkw = findTypeFromKeyword<BondedTypes>(header);
261             rtkw = findTypeFromKeyword<ReplaceType>(header);
262
263             if (!btkw.has_value() && !rtkw.has_value())
264             {
265                 tbptr->emplace_back(MoleculePatchDatabase());
266                 block = &tbptr->back();
267                 clearModificationBlock(block);
268                 block->name     = header;
269                 block->filebase = filebase;
270             }
271         }
272         else
273         {
274             if (block == nullptr)
275             {
276                 gmx_fatal(FARGS,
277                           "reading termini database: "
278                           "directive expected before line:\n%s\n"
279                           "This might be a file in an old format.",
280                           line);
281             }
282             /* this is not a header, so it must be data */
283             if (!btkw.has_value())
284             {
285                 /* this is a hack: add/rename/delete atoms */
286                 /* make space for hacks */
287                 block->hack.emplace_back(MoleculePatch());
288                 MoleculePatch* hack = &block->hack.back();
289
290                 /* get data */
291                 int n = 0;
292                 if (*rtkw == ReplaceType::Repl || *rtkw == ReplaceType::Del)
293                 {
294                     if (sscanf(line, "%s%n", buf, &n) != 1)
295                     {
296                         gmx_fatal(FARGS,
297                                   "Reading Termini Database '%s': "
298                                   "expected atom name on line\n%s",
299                                   fn,
300                                   line);
301                     }
302                     hack->oname = buf;
303                     /* we only replace or delete one atom at a time */
304                     hack->nr = 1;
305                 }
306                 else if (*rtkw == ReplaceType::Add)
307                 {
308                     read_ab(line, fn, hack);
309                     get_a_line(in, line, STRLEN);
310                 }
311                 else
312                 {
313                     gmx_fatal(FARGS,
314                               "unimplemented keyword number %d (%s:%d)",
315                               static_cast<int>(*rtkw),
316                               __FILE__,
317                               __LINE__);
318                 }
319                 if (*rtkw == ReplaceType::Repl || *rtkw == ReplaceType::Add)
320                 {
321                     hack->atom.emplace_back();
322                     read_atom(line + n, *rtkw == ReplaceType::Add, &hack->nname, &hack->atom.back(), atype, &hack->cgnr);
323                     if (hack->nname.empty())
324                     {
325                         if (!hack->oname.empty())
326                         {
327                             hack->nname = hack->oname;
328                         }
329                         else
330                         {
331                             gmx_fatal(FARGS,
332                                       "Reading Termini Database '%s': don't know which name the "
333                                       "new atom should have on line\n%s",
334                                       fn,
335                                       line);
336                         }
337                     }
338                 }
339             }
340             else if (*btkw >= BondedTypes::Bonds && *btkw < BondedTypes::Count)
341             {
342                 /* this is bonded data: bonds, angles, dihedrals or impropers */
343                 int n = 0;
344                 block->rb[*btkw].b.emplace_back();
345                 BondedInteraction* newBond = &block->rb[*btkw].b.back();
346                 for (int j = 0; j < enumValueToNumIAtoms(*btkw); j++)
347                 {
348                     int ni;
349                     if (sscanf(line + n, "%s%n", buf, &ni) == 1)
350                     {
351                         newBond->a[j] = buf;
352                     }
353                     else
354                     {
355                         gmx_fatal(FARGS,
356                                   "Reading Termini Database '%s': expected %d atom names (found "
357                                   "%d) on line\n%s",
358                                   fn,
359                                   enumValueToNumIAtoms(*btkw),
360                                   j - 1,
361                                   line);
362                     }
363                     n += ni;
364                 }
365                 strcpy(buf, "");
366                 sscanf(line + n, "%s", buf);
367                 newBond->s = buf;
368             }
369             else
370             {
371                 gmx_fatal(FARGS,
372                           "Reading Termini Database: Expecting a header at line\n"
373                           "%s",
374                           line);
375             }
376         }
377         get_a_line(in, line, STRLEN);
378     }
379
380     gmx_ffclose(in);
381 }
382
383 int read_ter_db(const char* ffdir, char ter, std::vector<MoleculePatchDatabase>* tbptr, PreprocessingAtomTypes* atype)
384 {
385     std::string ext = gmx::formatString(".%c.tdb", ter);
386
387     /* Search for termini database files.
388      * Do not generate an error when none are found.
389      */
390     std::vector<std::string> tdbf = fflib_search_file_end(ffdir, ext.c_str(), FALSE);
391     tbptr->clear();
392     for (const auto& filename : tdbf)
393     {
394         read_ter_db_file(filename.c_str(), tbptr, atype);
395     }
396
397     if (debug)
398     {
399         print_ter_db("new", ter, *tbptr, atype);
400     }
401
402     return tbptr->size();
403 }
404
405 std::vector<MoleculePatchDatabase*> filter_ter(gmx::ArrayRef<MoleculePatchDatabase> tb, const char* resname)
406 {
407     // TODO Four years later, no force fields have ever used this, so decide status of this feature
408     /* Since some force fields (e.g. OPLS) needs different
409      * atomtypes for different residues there could be a lot
410      * of entries in the databases for specific residues
411      * (e.g. GLY-NH3+, SER-NH3+, ALA-NH3+).
412      *
413      * To reduce the database size, we assume that a terminus specifier liker
414      *
415      * [ GLY|SER|ALA-NH3+ ]
416      *
417      * would cover all of the three residue types above.
418      * Wildcards (*,?) are not OK. Don't worry about e.g. GLU vs. GLUH since
419      * pdb2gmx only uses the first 3 letters when calling this routine.
420      *
421      * To automate this, this routines scans a list of termini
422      * for the residue name "resname" and returns an allocated list of
423      * pointers to the termini that could be applied to the
424      * residue in question. The variable pointed to by nret will
425      * contain the number of valid pointers in the list.
426      * Remember to free the list when you are done with it...
427      */
428
429     auto                                none_idx = tb.end();
430     std::vector<MoleculePatchDatabase*> list;
431
432     for (auto it = tb.begin(); it != tb.end(); it++)
433     {
434         const char* s     = it->name.c_str();
435         bool        found = false;
436         do
437         {
438             if (gmx::equalCaseInsensitive(resname, s, 3))
439             {
440                 found = true;
441                 list.push_back(&*it);
442             }
443             else
444             {
445                 /* advance to next |-separated field */
446                 s = strchr(s, '|');
447                 if (s != nullptr)
448                 {
449                     s++;
450                 }
451             }
452         } while (!found && s != nullptr);
453     }
454
455     /* All residue-specific termini have been added. We might have to fall
456      * back on generic termini, which are characterized by not having
457      * '-' in the name prior to the last position (which indicates charge).
458      * The [ None ] alternative is special since we don't want that
459      * to be the default, so we put it last in the list we return.
460      */
461     for (auto it = tb.begin(); it != tb.end(); it++)
462     {
463         const char* s = it->name.c_str();
464         if (gmx::equalCaseInsensitive("None", it->name))
465         {
466             none_idx = it;
467         }
468         else
469         {
470             /* Time to see if there's a generic terminus that matches.
471                Is there a hyphen? */
472             const char* c = strchr(s, '-');
473
474             /* A conjunction hyphen normally indicates a residue-specific
475                terminus, which is named like "GLY-COOH". A generic terminus
476                won't have a hyphen. */
477             bool bFoundAnyHyphen = (c != nullptr);
478             /* '-' as the last character indicates charge, so if that's
479                the only one found e.g. "COO-", then it was not a conjunction
480                hyphen, so this is a generic terminus */
481             bool bOnlyFoundChargeHyphen = (bFoundAnyHyphen && *(c + 1) == '\0');
482             /* Thus, "GLY-COO-" is not recognized as a generic terminus. */
483             bool bFoundGenericTerminus = !bFoundAnyHyphen || bOnlyFoundChargeHyphen;
484             if (bFoundGenericTerminus)
485             {
486                 /* Check that we haven't already added a residue-specific version
487                  * of this terminus.
488                  */
489                 auto found = std::find_if(list.begin(), list.end(), [&s](const MoleculePatchDatabase* b) {
490                     return strstr(b->name.c_str(), s) != nullptr;
491                 });
492                 if (found == list.end())
493                 {
494                     list.push_back(&*it);
495                 }
496             }
497         }
498     }
499     if (none_idx != tb.end())
500     {
501         list.push_back(&*none_idx);
502     }
503
504     return list;
505 }
506
507
508 MoleculePatchDatabase* choose_ter(gmx::ArrayRef<MoleculePatchDatabase*> tb, const char* title)
509 {
510     int sel, ret;
511
512     printf("%s\n", title);
513     int i = 0;
514     for (const auto& modification : tb)
515     {
516         bool bIsZwitterion = (0 == gmx_wcmatch("*ZWITTERION*", modification->name.c_str()));
517         printf("%2d: %s%s\n",
518                i,
519                modification->name.c_str(),
520                bIsZwitterion ? " (only use with zwitterions containing exactly one residue)" : "");
521         i++;
522     }
523     do
524     {
525         ret = fscanf(stdin, "%d", &sel);
526     } while ((ret != 1) || (sel < 0) || (sel >= tb.ssize()));
527
528     return tb[sel];
529 }