238aca95ca14c9def421e65878e08372c1316f96
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / readpull.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 <cassert>
40 #include <cstdlib>
41 #include <cstring>
42
43 #include "gromacs/domdec/localatomsetmanager.h"
44 #include "gromacs/fileio/readinp.h"
45 #include "gromacs/fileio/warninp.h"
46 #include "gromacs/gmxpreprocess/readir.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdlib/mdatoms.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/mdtypes/pull_params.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/pulling/pull.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
59
60
61 static void string2dvec(const char buf[], dvec nums)
62 {
63     double dum;
64
65     if (sscanf(buf, "%lf%lf%lf%lf", &nums[0], &nums[1], &nums[2], &dum) != 3)
66     {
67         gmx_fatal(FARGS, "Expected three numbers at input line %s", buf);
68     }
69 }
70
71 static void init_pull_group(t_pull_group* pg, const char* wbuf)
72 {
73     double d;
74     int    n;
75
76     pg->nweight = 0;
77     while (sscanf(wbuf, "%lf %n", &d, &n) == 1)
78     {
79         if (pg->nweight % 100 == 0)
80         {
81             srenew(pg->weight, pg->nweight + 100);
82         }
83         pg->weight[pg->nweight++] = d;
84         wbuf += n;
85     }
86 }
87
88 static void process_pull_dim(char* dim_buf, ivec dim, const t_pull_coord* pcrd)
89 {
90     int   ndim, d, nchar;
91     char *ptr, pulldim1[STRLEN];
92
93     ptr  = dim_buf;
94     ndim = 0;
95     for (d = 0; d < DIM; d++)
96     {
97         if (sscanf(ptr, "%s%n", pulldim1, &nchar) != 1)
98         {
99             gmx_fatal(FARGS, "Less than 3 pull dimensions given in pull_dim: '%s'", dim_buf);
100         }
101
102         if (gmx::equalCaseInsensitive(pulldim1, "N", 1))
103         {
104             dim[d] = 0;
105         }
106         else if (gmx::equalCaseInsensitive(pulldim1, "Y", 1))
107         {
108             dim[d] = 1;
109             ndim++;
110         }
111         else
112         {
113             gmx_fatal(FARGS, "Please use Y(ES) or N(O) for pull_dim only (not %s)", pulldim1);
114         }
115         ptr += nchar;
116     }
117     if (ndim == 0)
118     {
119         gmx_fatal(FARGS, "All entries in pull dim are N");
120     }
121     if ((pcrd->eGeom == epullgDIHEDRAL) && (ndim < 3))
122     {
123         gmx_fatal(FARGS, "Pull geometry dihedral is only useful with pull-dim = Y Y Y");
124     }
125     if ((pcrd->eGeom == epullgANGLE || pcrd->eGeom == epullgANGLEAXIS) && (ndim < 2))
126     {
127         gmx_fatal(FARGS,
128                   "Pull geometry %s is only useful with pull-dim = Y for at least 2 dimensions",
129                   EPULLGEOM(pcrd->eGeom));
130     }
131 }
132
133 static void init_pull_coord(t_pull_coord* pcrd,
134                             int           coord_index_for_output,
135                             char*         dim_buf,
136                             const char*   origin_buf,
137                             const char*   vec_buf,
138                             warninp_t     wi)
139 {
140     int  m;
141     dvec origin, vec;
142     char buf[STRLEN];
143
144     if (pcrd->eType == epullCONSTRAINT
145         && (pcrd->eGeom == epullgCYL || pcrd->eGeom == epullgDIRRELATIVE || pcrd->eGeom == epullgANGLE
146             || pcrd->eGeom == epullgANGLEAXIS || pcrd->eGeom == epullgDIHEDRAL))
147     {
148         gmx_fatal(FARGS,
149                   "Pulling of type %s can not be combined with geometry %s. Consider using pull "
150                   "type %s.",
151                   epull_names[pcrd->eType], epullg_names[pcrd->eGeom], epull_names[epullUMBRELLA]);
152     }
153
154     if (pcrd->eType == epullEXTERNAL)
155     {
156         if (pcrd->externalPotentialProvider[0] == '\0')
157         {
158             sprintf(buf,
159                     "The use of pull type '%s' for pull coordinate %d requires that the name of "
160                     "the module providing the potential external is set with the option %s%d%s",
161                     epull_names[pcrd->eType], coord_index_for_output, "pull-coord",
162                     coord_index_for_output, "-potential-provider");
163             warning_error(wi, buf);
164         }
165
166         if (pcrd->rate != 0)
167         {
168             sprintf(buf,
169                     "The use of pull type '%s' for pull coordinate %d requires that the pull rate "
170                     "is zero",
171                     epull_names[pcrd->eType], coord_index_for_output);
172             warning_error(wi, buf);
173         }
174
175         if (pcrd->eGeom == epullgCYL)
176         {
177             /* Warn the user of a PBC restriction, caused by the fact that
178              * there is no reference value with an external pull potential.
179              */
180             sprintf(buf,
181                     "With pull type '%s' and geometry '%s', the distance component along the "
182                     "cylinder axis between atoms in the cylinder group and the COM of the pull "
183                     "group should be smaller than half the box length",
184                     epull_names[pcrd->eType], epullg_names[pcrd->eGeom]);
185             warning_note(wi, buf);
186         }
187     }
188
189     process_pull_dim(dim_buf, pcrd->dim, pcrd);
190
191     string2dvec(origin_buf, origin);
192     if (pcrd->group[0] != 0 && dnorm(origin) > 0)
193     {
194         gmx_fatal(FARGS, "The pull origin can only be set with an absolute reference");
195     }
196
197     /* Check the given initial reference value and warn for dangerous values */
198     if (pcrd->eGeom == epullgDIST)
199     {
200         if (pcrd->bStart && pcrd->init < 0)
201         {
202             sprintf(buf,
203                     "The initial reference distance set by pull-coord-init is set to a negative "
204                     "value (%g) with geometry %s while distances need to be non-negative. "
205                     "This may work, since you have set pull-coord-start to 'yes' which modifies "
206                     "this value, but only for certain starting distances. "
207                     "If this is a mistake you may want to use geometry %s instead.",
208                     pcrd->init, EPULLGEOM(pcrd->eGeom), EPULLGEOM(epullgDIR));
209             warning(wi, buf);
210         }
211     }
212     else if (pcrd->eGeom == epullgANGLE || pcrd->eGeom == epullgANGLEAXIS)
213     {
214         if (pcrd->bStart && (pcrd->init < 0 || pcrd->init > 180))
215         {
216             /* This value of pcrd->init may be ok depending on pcrd->bStart which modifies pcrd->init later on */
217             sprintf(buf,
218                     "The initial reference angle set by pull-coord-init (%g) is outside of the "
219                     "allowed range [0, 180] degrees for geometry (%s). "
220                     "This may work, since you have set pull-coord-start to 'yes' which modifies "
221                     "this value, but only for certain starting angles.",
222                     pcrd->init, EPULLGEOM(pcrd->eGeom));
223             warning(wi, buf);
224         }
225     }
226     else if (pcrd->eGeom == epullgDIHEDRAL)
227     {
228         if (pcrd->bStart && (pcrd->init < -180 || pcrd->init > 180))
229         {
230             sprintf(buf,
231                     "The initial reference angle set by pull-coord-init (%g) is outside of the "
232                     "allowed range [-180, 180] degrees for geometry (%s). "
233                     "This may work, since you have set pull-coord-start to 'yes' which modifies "
234                     "this value, but only for certain starting angles.",
235                     pcrd->init, EPULLGEOM(pcrd->eGeom));
236             warning(wi, buf);
237         }
238     }
239
240     /* Check and set the pull vector */
241     clear_dvec(vec);
242     string2dvec(vec_buf, vec);
243
244     if (pcrd->eGeom == epullgDIR || pcrd->eGeom == epullgCYL || pcrd->eGeom == epullgDIRPBC
245         || pcrd->eGeom == epullgANGLEAXIS)
246     {
247         if (dnorm2(vec) == 0)
248         {
249             gmx_fatal(FARGS, "With pull geometry %s the pull vector can not be 0,0,0",
250                       epullg_names[pcrd->eGeom]);
251         }
252         for (int d = 0; d < DIM; d++)
253         {
254             if (vec[d] != 0 && pcrd->dim[d] == 0)
255             {
256                 gmx_fatal(FARGS,
257                           "pull-coord-vec has non-zero %c-component while pull_dim for the "
258                           "%c-dimension is set to N",
259                           'x' + d, 'x' + d);
260             }
261         }
262
263         /* Normalize the direction vector */
264         dsvmul(1 / dnorm(vec), vec, vec);
265     }
266     else /* This case is for are all the geometries where the pull vector is not used */
267     {
268         if (dnorm2(vec) > 0)
269         {
270             sprintf(buf,
271                     "A pull vector is given (%g  %g  %g) but will not be used with geometry %s. If "
272                     "you really want to use this "
273                     "vector, consider using geometry %s instead.",
274                     vec[0], vec[1], vec[2], EPULLGEOM(pcrd->eGeom),
275                     pcrd->eGeom == epullgANGLE ? EPULLGEOM(epullgANGLEAXIS) : EPULLGEOM(epullgDIR));
276             warning(wi, buf);
277         }
278     }
279     for (m = 0; m < DIM; m++)
280     {
281         pcrd->origin[m] = origin[m];
282         pcrd->vec[m]    = vec[m];
283     }
284 }
285
286 char** read_pullparams(std::vector<t_inpfile>* inp, pull_params_t* pull, warninp_t wi)
287 {
288     int    nscan, idum;
289     char** grpbuf;
290     char   buf[STRLEN];
291     char   provider[STRLEN], groups[STRLEN], dim_buf[STRLEN];
292     char   wbuf[STRLEN], origin_buf[STRLEN], vec_buf[STRLEN];
293
294     t_pull_group* pgrp;
295     t_pull_coord* pcrd;
296
297     /* read pull parameters */
298     printStringNoNewline(inp, "Cylinder radius for dynamic reaction force groups (nm)");
299     pull->cylinder_r     = get_ereal(inp, "pull-cylinder-r", 1.5, wi);
300     pull->constr_tol     = get_ereal(inp, "pull-constr-tol", 1E-6, wi);
301     pull->bPrintCOM      = (get_eeenum(inp, "pull-print-com", yesno_names, wi) != 0);
302     pull->bPrintRefValue = (get_eeenum(inp, "pull-print-ref-value", yesno_names, wi) != 0);
303     pull->bPrintComp     = (get_eeenum(inp, "pull-print-components", yesno_names, wi) != 0);
304     pull->nstxout        = get_eint(inp, "pull-nstxout", 50, wi);
305     pull->nstfout        = get_eint(inp, "pull-nstfout", 50, wi);
306     pull->bSetPbcRefToPrevStepCOM = (get_eeenum(inp, "pull-pbc-ref-prev-step-com", yesno_names, wi) != 0);
307     pull->bXOutAverage            = (get_eeenum(inp, "pull-xout-average", yesno_names, wi) != 0);
308     pull->bFOutAverage            = (get_eeenum(inp, "pull-fout-average", yesno_names, wi) != 0);
309     printStringNoNewline(inp, "Number of pull groups");
310     pull->ngroup = get_eint(inp, "pull-ngroups", 1, wi);
311     printStringNoNewline(inp, "Number of pull coordinates");
312     pull->ncoord = get_eint(inp, "pull-ncoords", 1, wi);
313
314     if (pull->ngroup < 1)
315     {
316         gmx_fatal(FARGS, "pull-ngroups should be >= 1");
317     }
318     /* We always add an absolute reference group (index 0), even if not used */
319     pull->ngroup += 1;
320
321     if (pull->ncoord < 1)
322     {
323         gmx_fatal(FARGS, "pull-ncoords should be >= 1");
324     }
325
326     snew(pull->group, pull->ngroup);
327
328     snew(pull->coord, pull->ncoord);
329
330     /* pull group options */
331     printStringNoNewline(inp, "Group and coordinate parameters");
332
333     /* Read the pull groups */
334     snew(grpbuf, pull->ngroup);
335     /* Group 0 is the absolute reference, we don't read anything for 0 */
336     for (int groupNum = 1; groupNum < pull->ngroup; groupNum++)
337     {
338         pgrp = &pull->group[groupNum];
339         snew(grpbuf[groupNum], STRLEN);
340         sprintf(buf, "pull-group%d-name", groupNum);
341         setStringEntry(inp, buf, grpbuf[groupNum], "");
342         sprintf(buf, "pull-group%d-weights", groupNum);
343         setStringEntry(inp, buf, wbuf, "");
344         sprintf(buf, "pull-group%d-pbcatom", groupNum);
345         pgrp->pbcatom = get_eint(inp, buf, 0, wi);
346
347         /* Initialize the pull group */
348         init_pull_group(pgrp, wbuf);
349     }
350
351     /* Read the pull coordinates */
352     for (int coordNum = 1; coordNum < pull->ncoord + 1; coordNum++)
353     {
354         pcrd = &pull->coord[coordNum - 1];
355         sprintf(buf, "pull-coord%d-type", coordNum);
356         pcrd->eType = get_eeenum(inp, buf, epull_names, wi);
357         sprintf(buf, "pull-coord%d-potential-provider", coordNum);
358         setStringEntry(inp, buf, provider, "");
359         pcrd->externalPotentialProvider = gmx_strdup(provider);
360         sprintf(buf, "pull-coord%d-geometry", coordNum);
361         pcrd->eGeom = get_eeenum(inp, buf, epullg_names, wi);
362         sprintf(buf, "pull-coord%d-groups", coordNum);
363         setStringEntry(inp, buf, groups, "");
364
365         switch (pcrd->eGeom)
366         {
367             case epullgDIHEDRAL: pcrd->ngroup = 6; break;
368             case epullgDIRRELATIVE:
369             case epullgANGLE: pcrd->ngroup = 4; break;
370             default: pcrd->ngroup = 2; break;
371         }
372
373         nscan = sscanf(groups, "%d %d %d %d %d %d %d", &pcrd->group[0], &pcrd->group[1],
374                        &pcrd->group[2], &pcrd->group[3], &pcrd->group[4], &pcrd->group[5], &idum);
375         if (nscan != pcrd->ngroup)
376         {
377             auto message =
378                     gmx::formatString("%s should contain %d pull group indices with geometry %s",
379                                       buf, pcrd->ngroup, epullg_names[pcrd->eGeom]);
380             set_warning_line(wi, nullptr, -1);
381             warning_error(wi, message);
382         }
383         for (int g = 0; g < pcrd->ngroup; g++)
384         {
385             if (pcrd->group[g] < 0 || pcrd->group[g] >= pull->ngroup)
386             {
387                 /* Quit with a fatal error to avoid invalid memory access */
388                 gmx_fatal(FARGS,
389                           "%s contains an invalid pull group %d, you should have %d <= group <= %d",
390                           buf, pcrd->group[g], 0, pull->ngroup - 1);
391             }
392         }
393
394         sprintf(buf, "pull-coord%d-dim", coordNum);
395         setStringEntry(inp, buf, dim_buf, "Y Y Y");
396         sprintf(buf, "pull-coord%d-origin", coordNum);
397         setStringEntry(inp, buf, origin_buf, "0.0 0.0 0.0");
398         sprintf(buf, "pull-coord%d-vec", coordNum);
399         setStringEntry(inp, buf, vec_buf, "0.0 0.0 0.0");
400         sprintf(buf, "pull-coord%d-start", coordNum);
401         pcrd->bStart = (get_eeenum(inp, buf, yesno_names, wi) != 0);
402         sprintf(buf, "pull-coord%d-init", coordNum);
403         pcrd->init = get_ereal(inp, buf, 0.0, wi);
404         sprintf(buf, "pull-coord%d-rate", coordNum);
405         pcrd->rate = get_ereal(inp, buf, 0.0, wi);
406         sprintf(buf, "pull-coord%d-k", coordNum);
407         pcrd->k = get_ereal(inp, buf, 0.0, wi);
408         sprintf(buf, "pull-coord%d-kB", coordNum);
409         pcrd->kB = get_ereal(inp, buf, pcrd->k, wi);
410
411         /* Initialize the pull coordinate */
412         init_pull_coord(pcrd, coordNum, dim_buf, origin_buf, vec_buf, wi);
413     }
414
415     return grpbuf;
416 }
417
418 void make_pull_groups(pull_params_t* pull, char** pgnames, const t_blocka* grps, char** gnames)
419 {
420     int           g, ig = -1, i;
421     t_pull_group* pgrp;
422
423     /* Absolute reference group (might not be used) is special */
424     pgrp                = &pull->group[0];
425     pgrp->nat           = 0;
426     pgrp->pbcatom       = -1;
427     pgrp->pbcatom_input = -1;
428
429     for (g = 1; g < pull->ngroup; g++)
430     {
431         pgrp                = &pull->group[g];
432         pgrp->pbcatom_input = pgrp->pbcatom;
433
434         if (strcmp(pgnames[g], "") == 0)
435         {
436             gmx_fatal(FARGS, "Pull option pull_group%d required by grompp has not been set.", g);
437         }
438
439         ig        = search_string(pgnames[g], grps->nr, gnames);
440         pgrp->nat = grps->index[ig + 1] - grps->index[ig];
441
442         fprintf(stderr, "Pull group %d '%s' has %d atoms\n", g, pgnames[g], pgrp->nat);
443
444         if (pgrp->nat == 0)
445         {
446             gmx_fatal(FARGS, "Pull group %d '%s' is empty", g, pgnames[g]);
447         }
448
449         snew(pgrp->ind, pgrp->nat);
450         for (i = 0; i < pgrp->nat; i++)
451         {
452             pgrp->ind[i] = grps->a[grps->index[ig] + i];
453         }
454
455         if (pgrp->nweight > 0 && pgrp->nweight != pgrp->nat)
456         {
457             gmx_fatal(FARGS,
458                       "Number of weights (%d) for pull group %d '%s' does not match the number of "
459                       "atoms (%d)",
460                       pgrp->nweight, g, pgnames[g], pgrp->nat);
461         }
462
463         if (pgrp->nat == 1)
464         {
465             /* No pbc is required for this group */
466             pgrp->pbcatom = -1;
467         }
468         else
469         {
470             if (pgrp->pbcatom > 0)
471             {
472                 pgrp->pbcatom -= 1;
473             }
474             else if (pgrp->pbcatom == 0)
475             {
476                 pgrp->pbcatom = pgrp->ind[(pgrp->nat - 1) / 2];
477             }
478             else
479             {
480                 /* Use cosine weighting */
481                 pgrp->pbcatom = -1;
482             }
483         }
484     }
485 }
486
487 void make_pull_coords(pull_params_t* pull)
488 {
489     int           c;
490     t_pull_coord* pcrd;
491
492     for (c = 0; c < pull->ncoord; c++)
493     {
494         pcrd = &pull->coord[c];
495
496         if (pcrd->group[0] < 0 || pcrd->group[0] >= pull->ngroup || pcrd->group[1] < 0
497             || pcrd->group[1] >= pull->ngroup)
498         {
499             gmx_fatal(FARGS,
500                       "Pull group index in pull-coord%d-groups out of range, should be between %d "
501                       "and %d",
502                       c + 1, 0, pull->ngroup + 1);
503         }
504
505         if (pcrd->group[0] == pcrd->group[1])
506         {
507             gmx_fatal(FARGS, "Identical pull group indices in pull-coord%d-groups", c + 1);
508         }
509
510         if (pcrd->eGeom == epullgCYL)
511         {
512             if (pull->group[pcrd->group[0]].nweight > 0)
513             {
514                 gmx_fatal(
515                         FARGS,
516                         "Weights are not supported for the reference group with cylinder pulling");
517             }
518         }
519     }
520 }
521
522 pull_t* set_pull_init(t_inputrec* ir, const gmx_mtop_t* mtop, rvec* x, matrix box, real lambda, warninp_t wi)
523 {
524     pull_params_t* pull;
525     pull_t*        pull_work;
526     t_pbc          pbc;
527     int            c;
528     double         t_start;
529
530     pull = ir->pull;
531     gmx::LocalAtomSetManager atomSets;
532     pull_work    = init_pull(nullptr, pull, ir, mtop, nullptr, &atomSets, lambda);
533     auto mdAtoms = gmx::makeMDAtoms(nullptr, *mtop, *ir, false);
534     auto md      = mdAtoms->mdatoms();
535     atoms2md(mtop, ir, -1, nullptr, mtop->natoms, mdAtoms.get());
536     if (ir->efep)
537     {
538         update_mdatoms(md, lambda);
539     }
540
541     set_pbc(&pbc, ir->ePBC, box);
542
543     t_start = ir->init_t + ir->init_step * ir->delta_t;
544
545     if (pull->bSetPbcRefToPrevStepCOM)
546     {
547         initPullComFromPrevStep(nullptr, pull_work, md, &pbc, x);
548     }
549     pull_calc_coms(nullptr, pull_work, md, &pbc, t_start, x, nullptr);
550
551     for (int g = 0; g < pull->ngroup; g++)
552     {
553         bool groupObeysPbc = pullCheckPbcWithinGroup(
554                 *pull_work, gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec*>(x), mtop->natoms),
555                 pbc, g, c_pullGroupSmallGroupThreshold);
556         if (!groupObeysPbc)
557         {
558             char buf[STRLEN];
559             if (pull->group[g].pbcatom_input == 0)
560             {
561                 sprintf(buf,
562                         "When the maximum distance from a pull group reference atom to other atoms "
563                         "in the "
564                         "group is larger than %g times half the box size a centrally placed "
565                         "atom should be chosen as pbcatom. Pull group %d is larger than that and "
566                         "does not have "
567                         "a specific atom selected as reference atom.",
568                         c_pullGroupSmallGroupThreshold, g);
569                 warning_error(wi, buf);
570             }
571             else if (!pull->bSetPbcRefToPrevStepCOM)
572             {
573                 sprintf(buf,
574                         "The maximum distance from the chosen PBC atom (%d) of pull group %d to "
575                         "other "
576                         "atoms in the group is larger than %g times half the box size. "
577                         "Set the pull-pbc-ref-prev-step-com option to yes.",
578                         pull->group[g].pbcatom + 1, g, c_pullGroupSmallGroupThreshold);
579                 warning_error(wi, buf);
580             }
581         }
582         if (groupObeysPbc)
583         {
584             groupObeysPbc = pullCheckPbcWithinGroup(
585                     *pull_work, gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec*>(x), mtop->natoms),
586                     pbc, g, c_pullGroupPbcMargin);
587             if (!groupObeysPbc)
588             {
589                 char buf[STRLEN];
590                 sprintf(buf,
591                         "Pull group %d has atoms at a distance larger than %g times half the box "
592                         "size from the PBC atom (%d). "
593                         "If atoms are or will more beyond half the box size from the PBC atom, the "
594                         "COM will be ill defined.",
595                         g, c_pullGroupPbcMargin, pull->group[g].pbcatom + 1);
596                 set_warning_line(wi, nullptr, -1);
597                 warning(wi, buf);
598             }
599         }
600     }
601
602     fprintf(stderr, "Pull group  natoms  pbc atom  distance at start  reference at t=0\n");
603     for (c = 0; c < pull->ncoord; c++)
604     {
605         t_pull_coord* pcrd;
606         t_pull_group *pgrp0, *pgrp1;
607         double        value;
608         real          init = 0;
609
610         pcrd = &pull->coord[c];
611
612         pgrp0 = &pull->group[pcrd->group[0]];
613         pgrp1 = &pull->group[pcrd->group[1]];
614         fprintf(stderr, "%8d  %8d  %8d\n", pcrd->group[0], pgrp0->nat, pgrp0->pbcatom + 1);
615         fprintf(stderr, "%8d  %8d  %8d ", pcrd->group[1], pgrp1->nat, pgrp1->pbcatom + 1);
616
617         if (pcrd->bStart)
618         {
619             init       = pcrd->init;
620             pcrd->init = 0;
621         }
622
623         value = get_pull_coord_value(pull_work, c, &pbc);
624
625         value *= pull_conversion_factor_internal2userinput(pcrd);
626         fprintf(stderr, " %10.3f %s", value, pull_coordinate_units(pcrd));
627
628         if (pcrd->bStart)
629         {
630             pcrd->init = value + init;
631         }
632
633         if (pcrd->eGeom == epullgDIST)
634         {
635             if (pcrd->init < 0)
636             {
637                 gmx_fatal(FARGS,
638                           "The initial pull distance (%g) needs to be non-negative with geometry "
639                           "%s. If you want a signed distance, use geometry %s instead.",
640                           pcrd->init, EPULLGEOM(pcrd->eGeom), EPULLGEOM(epullgDIR));
641             }
642
643             /* TODO: With a positive init but a negative rate things could still
644              * go wrong, but it might be fine if you don't pull too far.
645              * We should give a warning or note when there is only one pull dim
646              * active, since that is usually the problematic case when you should
647              * be using direction. We will do this later, since an already planned
648              * generalization of the pull code makes pull dim available here.
649              */
650         }
651         else if (pcrd->eGeom == epullgANGLE || pcrd->eGeom == epullgANGLEAXIS)
652         {
653             if (pcrd->init < 0 || pcrd->init > 180)
654             {
655                 gmx_fatal(FARGS,
656                           "The initial pull reference angle (%g) is outside of the allowed range "
657                           "[0, 180] degrees.",
658                           pcrd->init);
659             }
660         }
661         else if (pcrd->eGeom == epullgDIHEDRAL)
662         {
663             if (pcrd->init < -180 || pcrd->init > 180)
664             {
665                 gmx_fatal(FARGS,
666                           "The initial pull reference angle (%g) is outside of the allowed range "
667                           "[-180, 180] degrees.",
668                           pcrd->init);
669             }
670         }
671
672
673         fprintf(stderr, "     %10.3f %s\n", pcrd->init, pull_coordinate_units(pcrd));
674     }
675
676     return pull_work;
677 }