Do not include headers related to ObservablesHistory
[alexxy/gromacs.git] / src / programs / mdrun / runner.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) 2011,2012,2013,2014,2015,2016,2017, 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 /*! \internal \file
38  *
39  * \brief Implements the MD runner routine calling all integrators.
40  *
41  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42  * \ingroup module_mdlib
43  */
44 #include "gmxpre.h"
45
46 #include "runner.h"
47
48 #include "config.h"
49
50 #include <cassert>
51 #include <csignal>
52 #include <cstdlib>
53 #include <cstring>
54
55 #include <algorithm>
56
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/ewald/ewald-utils.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/checkpoint.h"
63 #include "gromacs/fileio/oenv.h"
64 #include "gromacs/fileio/tpxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gpu_utils/gpu_utils.h"
67 #include "gromacs/hardware/cpuinfo.h"
68 #include "gromacs/hardware/detecthardware.h"
69 #include "gromacs/hardware/printhardware.h"
70 #include "gromacs/listed-forces/disre.h"
71 #include "gromacs/listed-forces/orires.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/utilities.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/mdlib/calc_verletbuf.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/force.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/gmx_omp_nthreads.h"
80 #include "gromacs/mdlib/integrator.h"
81 #include "gromacs/mdlib/main.h"
82 #include "gromacs/mdlib/md_support.h"
83 #include "gromacs/mdlib/mdatoms.h"
84 #include "gromacs/mdlib/mdrun.h"
85 #include "gromacs/mdlib/minimize.h"
86 #include "gromacs/mdlib/nb_verlet.h"
87 #include "gromacs/mdlib/nbnxn_search.h"
88 #include "gromacs/mdlib/nbnxn_tuning.h"
89 #include "gromacs/mdlib/qmmm.h"
90 #include "gromacs/mdlib/sighandler.h"
91 #include "gromacs/mdlib/sim_util.h"
92 #include "gromacs/mdlib/tpi.h"
93 #include "gromacs/mdrunutility/mdmodules.h"
94 #include "gromacs/mdrunutility/threadaffinity.h"
95 #include "gromacs/mdtypes/commrec.h"
96 #include "gromacs/mdtypes/inputrec.h"
97 #include "gromacs/mdtypes/md_enums.h"
98 #include "gromacs/mdtypes/observableshistory.h"
99 #include "gromacs/mdtypes/state.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/pulling/pull.h"
102 #include "gromacs/pulling/pull_rotation.h"
103 #include "gromacs/taskassignment/hardwareassign.h"
104 #include "gromacs/taskassignment/resourcedivision.h"
105 #include "gromacs/timing/wallcycle.h"
106 #include "gromacs/topology/mtop_util.h"
107 #include "gromacs/trajectory/trajectoryframe.h"
108 #include "gromacs/utility/cstringutil.h"
109 #include "gromacs/utility/exceptions.h"
110 #include "gromacs/utility/fatalerror.h"
111 #include "gromacs/utility/filestream.h"
112 #include "gromacs/utility/gmxassert.h"
113 #include "gromacs/utility/gmxmpi.h"
114 #include "gromacs/utility/logger.h"
115 #include "gromacs/utility/loggerbuilder.h"
116 #include "gromacs/utility/pleasecite.h"
117 #include "gromacs/utility/programcontext.h"
118 #include "gromacs/utility/smalloc.h"
119
120 #include "deform.h"
121 #include "md.h"
122 #include "membed.h"
123 #include "repl_ex.h"
124
125 #ifdef GMX_FAHCORE
126 #include "corewrap.h"
127 #endif
128
129 //! First step used in pressure scaling
130 gmx_int64_t         deform_init_init_step_tpx;
131 //! Initial box for pressure scaling
132 matrix              deform_init_box_tpx;
133 //! MPI variable for use in pressure scaling
134 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
135
136 #if GMX_THREAD_MPI
137 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
138  * the number of threads will get lowered.
139  */
140 #define MIN_ATOMS_PER_MPI_THREAD    90
141 #define MIN_ATOMS_PER_GPU           900
142
143 namespace gmx
144 {
145
146 void Mdrunner::reinitializeOnSpawnedThread()
147 {
148     // TODO This duplication is formally necessary if any thread might
149     // modify any memory in fnm or the pointers it contains. If the
150     // contents are ever provably const, then we can remove this
151     // allocation (and memory leak).
152     // TODO This should probably become part of a copy constructor for
153     // Mdrunner.
154     fnm = dup_tfn(nfile, fnm);
155
156     cr  = reinitialize_commrec_for_this_thread(cr);
157
158     if (!MASTER(cr))
159     {
160         // Only the master rank writes to the log files
161         fplog = nullptr;
162     }
163 }
164
165 /*! \brief The callback used for running on spawned threads.
166  *
167  * Obtains the pointer to the master mdrunner object from the one
168  * argument permitted to the thread-launch API call, copies it to make
169  * a new runner for this thread, reinitializes necessary data, and
170  * proceeds to the simulation. */
171 static void mdrunner_start_fn(void *arg)
172 {
173     try
174     {
175         auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
176         /* copy the arg list to make sure that it's thread-local. This
177            doesn't copy pointed-to items, of course, but those are all
178            const. */
179         gmx::Mdrunner mdrunner = *masterMdrunner;
180         mdrunner.reinitializeOnSpawnedThread();
181         mdrunner.mdrunner();
182     }
183     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
184 }
185
186
187 /*! \brief Start thread-MPI threads.
188  *
189  * Called by mdrunner() to start a specific number of threads
190  * (including the main thread) for thread-parallel runs. This in turn
191  * calls mdrunner() for each thread. All options are the same as for
192  * mdrunner(). */
193 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch)
194 {
195
196     /* first check whether we even need to start tMPI */
197     if (numThreadsToLaunch < 2)
198     {
199         return cr;
200     }
201
202     gmx::Mdrunner spawnedMdrunner = *this;
203     // TODO This duplication is formally necessary if any thread might
204     // modify any memory in fnm or the pointers it contains. If the
205     // contents are ever provably const, then we can remove this
206     // allocation (and memory leak).
207     // TODO This should probably become part of a copy constructor for
208     // Mdrunner.
209     spawnedMdrunner.fnm = dup_tfn(this->nfile, fnm);
210
211     /* now spawn new threads that start mdrunner_start_fn(), while
212        the main thread returns, we set thread affinity later */
213     if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
214                      mdrunner_start_fn, static_cast<void*>(&spawnedMdrunner)) != TMPI_SUCCESS)
215     {
216         GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
217     }
218
219     return reinitialize_commrec_for_this_thread(cr);
220 }
221
222 }      // namespace
223
224 #endif /* GMX_THREAD_MPI */
225
226 /*! \brief Initialize variables for Verlet scheme simulation */
227 static void prepare_verlet_scheme(FILE                           *fplog,
228                                   t_commrec                      *cr,
229                                   t_inputrec                     *ir,
230                                   int                             nstlist_cmdline,
231                                   const gmx_mtop_t               *mtop,
232                                   const matrix                    box,
233                                   bool                            makeGpuPairList,
234                                   const gmx::CpuInfo             &cpuinfo)
235 {
236     /* For NVE simulations, we will retain the initial list buffer */
237     if (EI_DYNAMICS(ir->eI) &&
238         ir->verletbuf_tol > 0 &&
239         !(EI_MD(ir->eI) && ir->etc == etcNO))
240     {
241         /* Update the Verlet buffer size for the current run setup */
242         verletbuf_list_setup_t ls;
243         real                   rlist_new;
244
245         /* Here we assume SIMD-enabled kernels are being used. But as currently
246          * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
247          * and 4x2 gives a larger buffer than 4x4, this is ok.
248          */
249         verletbuf_get_list_setup(true, makeGpuPairList, &ls);
250
251         calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1, -1, &ls, nullptr, &rlist_new);
252
253         if (rlist_new != ir->rlist)
254         {
255             if (fplog != nullptr)
256             {
257                 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
258                         ir->rlist, rlist_new,
259                         ls.cluster_size_i, ls.cluster_size_j);
260             }
261             ir->rlist     = rlist_new;
262         }
263     }
264
265     if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
266     {
267         gmx_fatal(FARGS, "Can not set nstlist without %s",
268                   !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
269     }
270
271     if (EI_DYNAMICS(ir->eI))
272     {
273         /* Set or try nstlist values */
274         increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
275     }
276 }
277
278 /*! \brief Override the nslist value in inputrec
279  *
280  * with value passed on the command line (if any)
281  */
282 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
283                                     gmx_int64_t          nsteps_cmdline,
284                                     t_inputrec          *ir)
285 {
286     assert(ir);
287
288     /* override with anything else than the default -2 */
289     if (nsteps_cmdline > -2)
290     {
291         char sbuf_steps[STEPSTRSIZE];
292         char sbuf_msg[STRLEN];
293
294         ir->nsteps = nsteps_cmdline;
295         if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
296         {
297             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
298                     gmx_step_str(nsteps_cmdline, sbuf_steps),
299                     fabs(nsteps_cmdline*ir->delta_t));
300         }
301         else
302         {
303             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
304                     gmx_step_str(nsteps_cmdline, sbuf_steps));
305         }
306
307         GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
308     }
309     else if (nsteps_cmdline < -2)
310     {
311         gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
312                   nsteps_cmdline);
313     }
314     /* Do nothing if nsteps_cmdline == -2 */
315 }
316
317 namespace gmx
318 {
319
320 //! Halt the run if there are inconsistences between user choices to run with GPUs and/or hardware detection.
321 static void exitIfCannotForceGpuRun(bool                requirePhysicalGpu,
322                                     EmulateGpuNonbonded emulateGpuNonbonded,
323                                     bool                useVerletScheme,
324                                     bool                compatibleGpusFound)
325 {
326     /* Was GPU acceleration either explicitly (-nb gpu) or implicitly
327      * (gpu ID passed) requested? */
328     if (!requirePhysicalGpu)
329     {
330         return;
331     }
332
333     if (GMX_GPU == GMX_GPU_NONE)
334     {
335         gmx_fatal(FARGS, "GPU acceleration requested, but %s was compiled without GPU support!",
336                   gmx::getProgramContext().displayName());
337     }
338
339     if (emulateGpuNonbonded == EmulateGpuNonbonded::Yes)
340     {
341         gmx_fatal(FARGS, "GPU emulation cannot be requested together with GPU acceleration!");
342     }
343
344     if (!useVerletScheme)
345     {
346         gmx_fatal(FARGS, "GPU acceleration requested, but can't be used without cutoff-scheme=Verlet");
347     }
348
349     if (!compatibleGpusFound)
350     {
351         gmx_fatal(FARGS, "GPU acceleration requested, but no compatible GPUs were detected.");
352     }
353 }
354
355 /*! \brief Return whether GPU acceleration is useful with the given settings.
356  *
357  * If not, logs a message about falling back to CPU code. */
358 static bool gpuAccelerationIsUseful(const MDLogger   &mdlog,
359                                     const t_inputrec *ir,
360                                     bool              doRerun)
361 {
362     if (doRerun && ir->opts.ngener > 1)
363     {
364         /* Rerun execution time is dominated by I/O and pair search,
365          * so GPUs are not very useful, plus they do not support more
366          * than one energy group. If the user requested GPUs
367          * explicitly, a fatal error is given later.  With non-reruns,
368          * we fall back to a single whole-of system energy group
369          * (which runs much faster than a multiple-energy-groups
370          * implementation would), and issue a note in the .log
371          * file. Users can re-run if they want the information. */
372         GMX_LOG(mdlog.warning).asParagraph().appendText("Multiple energy groups is not implemented for GPUs, so is not useful for this rerun, so falling back to the CPU");
373         return false;
374     }
375
376     return true;
377 }
378
379 //! \brief Return the correct integrator function.
380 static integrator_t *my_integrator(unsigned int ei)
381 {
382     switch (ei)
383     {
384         case eiMD:
385         case eiBD:
386         case eiSD1:
387         case eiVV:
388         case eiVVAK:
389             if (!EI_DYNAMICS(ei))
390             {
391                 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
392             }
393             return do_md;
394         case eiSteep:
395             return do_steep;
396         case eiCG:
397             return do_cg;
398         case eiNM:
399             return do_nm;
400         case eiLBFGS:
401             return do_lbfgs;
402         case eiTPI:
403         case eiTPIC:
404             if (!EI_TPI(ei))
405             {
406                 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
407             }
408             return do_tpi;
409         case eiSD2_REMOVED:
410             GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
411         default:
412             GMX_THROW(APIError("Non existing integrator selected"));
413     }
414 }
415
416 //! Initializes the logger for mdrun.
417 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
418 {
419     gmx::LoggerBuilder builder;
420     if (fplog != nullptr)
421     {
422         builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
423     }
424     if (cr == nullptr || SIMMASTER(cr))
425     {
426         builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
427                                 &gmx::TextOutputFile::standardError());
428     }
429     return builder.build();
430 }
431
432 int Mdrunner::mdrunner()
433 {
434     matrix                    box;
435     gmx_ddbox_t               ddbox = {0};
436     int                       npme_major, npme_minor;
437     t_nrnb                   *nrnb;
438     gmx_mtop_t               *mtop          = nullptr;
439     t_mdatoms                *mdatoms       = nullptr;
440     t_forcerec               *fr            = nullptr;
441     t_fcdata                 *fcd           = nullptr;
442     real                      ewaldcoeff_q  = 0;
443     real                      ewaldcoeff_lj = 0;
444     struct gmx_pme_t        **pmedata       = nullptr;
445     gmx_vsite_t              *vsite         = nullptr;
446     gmx_constr_t              constr;
447     int                       nChargePerturbed = -1, nTypePerturbed = 0, status;
448     gmx_wallcycle_t           wcycle;
449     gmx_walltime_accounting_t walltime_accounting = nullptr;
450     int                       rc;
451     gmx_int64_t               reset_counters;
452     int                       nthreads_pme = 1;
453     gmx_membed_t *            membed       = nullptr;
454     gmx_hw_info_t            *hwinfo       = nullptr;
455
456     /* CAUTION: threads may be started later on in this function, so
457        cr doesn't reflect the final parallel state right now */
458     gmx::MDModules mdModules;
459     t_inputrec     inputrecInstance;
460     t_inputrec    *inputrec = &inputrecInstance;
461     snew(mtop, 1);
462
463     if (mdrunOptions.continuationOptions.appendFiles)
464     {
465         fplog = nullptr;
466     }
467
468     bool doMembed = opt2bSet("-membed", nfile, fnm);
469     bool doRerun  = mdrunOptions.rerun;
470
471     /* Handle GPU-related user options. Later, we check consistency
472      * with things like whether support is compiled, or tMPI thread
473      * count. */
474     EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
475                                                EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
476     bool                forceUseCpu           = (strncmp(nbpu_opt, "cpu", 3) == 0);
477     if (!hw_opt.gpuIdTaskAssignment.empty() && forceUseCpu)
478     {
479         gmx_fatal(FARGS, "GPU IDs were specified, and short-ranged interactions were assigned to the CPU. Make no more than one of these choices.");
480     }
481     bool forceUsePhysicalGpu = (strncmp(nbpu_opt, "gpu", 3) == 0) || !hw_opt.gpuIdTaskAssignment.empty();
482     bool tryUsePhysicalGpu   = (strncmp(nbpu_opt, "auto", 4) == 0) && hw_opt.gpuIdTaskAssignment.empty() && (emulateGpuNonbonded == EmulateGpuNonbonded::No);
483     GMX_RELEASE_ASSERT(!(forceUsePhysicalGpu && tryUsePhysicalGpu), "Must either force use of "
484                        "GPUs for short-ranged interactions, or try to use them, not both.");
485
486     // Here we assume that SIMMASTER(cr) does not change even after the
487     // threads are started.
488     gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
489     gmx::MDLogger    mdlog(logOwner.logger());
490
491     hwinfo = gmx_detect_hardware(mdlog, cr);
492
493     gmx_print_detected_hardware(fplog, cr, mdlog, hwinfo);
494
495     if (fplog != nullptr)
496     {
497         /* Print references after all software/hardware printing */
498         please_cite(fplog, "Abraham2015");
499         please_cite(fplog, "Pall2015");
500         please_cite(fplog, "Pronk2013");
501         please_cite(fplog, "Hess2008b");
502         please_cite(fplog, "Spoel2005a");
503         please_cite(fplog, "Lindahl2001a");
504         please_cite(fplog, "Berendsen95a");
505     }
506
507     std::unique_ptr<t_state> globalState;
508
509     if (SIMMASTER(cr))
510     {
511         /* Only the master rank has the global state */
512         globalState = std::unique_ptr<t_state>(new t_state);
513
514         /* Read (nearly) all data required for the simulation */
515         read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, globalState.get(), mtop);
516
517         exitIfCannotForceGpuRun(forceUsePhysicalGpu,
518                                 emulateGpuNonbonded,
519                                 inputrec->cutoff_scheme == ecutsVERLET,
520                                 compatibleGpusFound(hwinfo->gpu_info));
521
522         if (inputrec->cutoff_scheme == ecutsVERLET)
523         {
524             /* TODO This logic could run later, e.g. before -npme -1
525                is handled. If inputrec has already been communicated,
526                then the resulting tryUsePhysicalGpu does not need to
527                be communicated. */
528             if ((tryUsePhysicalGpu || forceUsePhysicalGpu) &&
529                 !gpuAccelerationIsUseful(mdlog, inputrec, doRerun))
530             {
531                 /* Fallback message printed by nbnxn_acceleration_supported */
532                 if (forceUsePhysicalGpu)
533                 {
534                     gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
535                 }
536                 tryUsePhysicalGpu = false;
537             }
538         }
539         else
540         {
541             if (nstlist_cmdline > 0)
542             {
543                 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
544             }
545
546             if (compatibleGpusFound(hwinfo->gpu_info))
547             {
548                 GMX_LOG(mdlog.warning).asParagraph().appendText(
549                         "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
550                         "      To use a GPU, set the mdp option: cutoff-scheme = Verlet");
551             }
552             tryUsePhysicalGpu = false;
553         }
554     }
555     bool nonbondedOnGpu = (tryUsePhysicalGpu || forceUsePhysicalGpu) && compatibleGpusFound(hwinfo->gpu_info);
556
557     /* Check and update the hardware options for internal consistency */
558     check_and_update_hw_opt_1(&hw_opt, cr, domdecOptions.numPmeRanks);
559
560     /* Early check for externally set process affinity. */
561     gmx_check_thread_affinity_set(mdlog, cr,
562                                   &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
563
564 #if GMX_THREAD_MPI
565     if (SIMMASTER(cr))
566     {
567         if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
568         {
569             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
570         }
571
572         /* Since the master knows the cut-off scheme, update hw_opt for this.
573          * This is done later for normal MPI and also once more with tMPI
574          * for all tMPI ranks.
575          */
576         check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
577
578         /* Determine how many thread-MPI ranks to start.
579          *
580          * TODO Over-writing the user-supplied value here does
581          * prevent any possible subsequent checks from working
582          * correctly. */
583         hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
584                                                 &hw_opt,
585                                                 domdecOptions.numPmeRanks,
586                                                 nonbondedOnGpu,
587                                                 inputrec, mtop,
588                                                 mdlog,
589                                                 doMembed);
590
591         // Now start the threads for thread MPI.
592         cr = spawnThreads(hw_opt.nthreads_tmpi);
593         /* The main thread continues here with a new cr. We don't deallocate
594            the old cr because other threads may still be reading it. */
595         // TODO Both master and spawned threads call dup_tfn and
596         // reinitialize_commrec_for_this_thread. Find a way to express
597         // this better.
598     }
599 #endif
600     /* END OF CAUTION: cr is now reliable */
601
602     if (PAR(cr))
603     {
604         /* now broadcast everything to the non-master nodes/threads: */
605         init_parallel(cr, inputrec, mtop);
606
607         gmx_bcast_sim(sizeof(nonbondedOnGpu), &nonbondedOnGpu, cr);
608     }
609     // TODO: Error handling
610     mdModules.assignOptionsToModules(*inputrec->params, nullptr);
611
612     if (fplog != nullptr)
613     {
614         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
615         fprintf(fplog, "\n");
616     }
617
618     if (SIMMASTER(cr))
619     {
620         /* now make sure the state is initialized and propagated */
621         set_state_entries(globalState.get(), inputrec);
622     }
623
624     /* NM and TPI parallelize over force/energy calculations, not atoms,
625      * so we need to initialize and broadcast the global state.
626      */
627     if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
628     {
629         if (!MASTER(cr))
630         {
631             globalState = std::unique_ptr<t_state>(new t_state);
632         }
633         broadcastStateWithoutDynamics(cr, globalState.get());
634     }
635
636     /* A parallel command line option consistency check that we can
637        only do after any threads have started. */
638     if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
639                      domdecOptions.numCells[YY] > 1 ||
640                      domdecOptions.numCells[ZZ] > 1 ||
641                      domdecOptions.numPmeRanks > 0))
642     {
643         gmx_fatal(FARGS,
644                   "The -dd or -npme option request a parallel simulation, "
645 #if !GMX_MPI
646                   "but %s was compiled without threads or MPI enabled"
647 #else
648 #if GMX_THREAD_MPI
649                   "but the number of MPI-threads (option -ntmpi) is not set or is 1"
650 #else
651                   "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
652 #endif
653 #endif
654                   , output_env_get_program_display_name(oenv)
655                   );
656     }
657
658     if (doRerun &&
659         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
660     {
661         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
662     }
663
664     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
665     {
666         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
667     }
668
669     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
670     {
671         if (domdecOptions.numPmeRanks > 0)
672         {
673             gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
674                                  "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
675         }
676
677         domdecOptions.numPmeRanks = 0;
678     }
679
680     if (nonbondedOnGpu && domdecOptions.numPmeRanks < 0)
681     {
682         /* With GPUs we don't automatically use PME-only ranks. PME ranks can
683          * improve performance with many threads per GPU, since our OpenMP
684          * scaling is bad, but it's difficult to automate the setup.
685          */
686         domdecOptions.numPmeRanks = 0;
687     }
688
689 #ifdef GMX_FAHCORE
690     if (MASTER(cr))
691     {
692         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
693     }
694 #endif
695
696     /* NMR restraints must be initialized before load_checkpoint,
697      * since with time averaging the history is added to t_state.
698      * For proper consistency check we therefore need to extend
699      * t_state here.
700      * So the PME-only nodes (if present) will also initialize
701      * the distance restraints.
702      */
703     snew(fcd, 1);
704
705     /* This needs to be called before read_checkpoint to extend the state */
706     init_disres(fplog, mtop, inputrec, cr, fcd, globalState.get(), replExParams.exchangeInterval > 0);
707
708     init_orires(fplog, mtop, inputrec, cr, globalState.get(), &(fcd->orires));
709
710     if (inputrecDeform(inputrec))
711     {
712         /* Store the deform reference box before reading the checkpoint */
713         if (SIMMASTER(cr))
714         {
715             copy_mat(globalState->box, box);
716         }
717         if (PAR(cr))
718         {
719             gmx_bcast(sizeof(box), box, cr);
720         }
721         /* Because we do not have the update struct available yet
722          * in which the reference values should be stored,
723          * we store them temporarily in static variables.
724          * This should be thread safe, since they are only written once
725          * and with identical values.
726          */
727         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
728         deform_init_init_step_tpx = inputrec->init_step;
729         copy_mat(box, deform_init_box_tpx);
730         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
731     }
732
733     ObservablesHistory   observablesHistory = {};
734
735     ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
736
737     if (continuationOptions.startedFromCheckpoint)
738     {
739         /* Check if checkpoint file exists before doing continuation.
740          * This way we can use identical input options for the first and subsequent runs...
741          */
742         gmx_bool bReadEkin;
743
744         load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
745                         cr, domdecOptions.numCells,
746                         inputrec, globalState.get(),
747                         &bReadEkin, &observablesHistory,
748                         continuationOptions.appendFiles,
749                         continuationOptions.appendFilesOptionSet,
750                         mdrunOptions.reproducible);
751
752         if (bReadEkin)
753         {
754             continuationOptions.haveReadEkin = true;
755         }
756     }
757
758     if (SIMMASTER(cr) && continuationOptions.appendFiles)
759     {
760         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
761                      continuationOptions.appendFiles, &fplog);
762         logOwner = buildLogger(fplog, nullptr);
763         mdlog    = logOwner.logger();
764     }
765
766     /* override nsteps with value set on the commamdline */
767     override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
768
769     if (SIMMASTER(cr))
770     {
771         copy_mat(globalState->box, box);
772     }
773
774     if (PAR(cr))
775     {
776         gmx_bcast(sizeof(box), box, cr);
777     }
778
779     /* Update rlist and nstlist. */
780     if (inputrec->cutoff_scheme == ecutsVERLET)
781     {
782         prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, mtop, box,
783                               nonbondedOnGpu || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
784     }
785
786     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
787                      inputrec->eI == eiNM))
788     {
789         const rvec *xOnMaster = (SIMMASTER(cr) ? as_rvec_array(globalState->x.data()) : nullptr);
790
791         cr->dd = init_domain_decomposition(fplog, cr, domdecOptions, mdrunOptions,
792                                            mtop, inputrec,
793                                            box, xOnMaster,
794                                            &ddbox, &npme_major, &npme_minor);
795         // Note that local state still does not exist yet.
796     }
797     else
798     {
799         /* PME, if used, is done on all nodes with 1D decomposition */
800         cr->npmenodes = 0;
801         cr->duty      = (DUTY_PP | DUTY_PME);
802         npme_major    = 1;
803         npme_minor    = 1;
804
805         if (inputrec->ePBC == epbcSCREW)
806         {
807             gmx_fatal(FARGS,
808                       "pbc=%s is only implemented with domain decomposition",
809                       epbc_names[inputrec->ePBC]);
810         }
811     }
812
813     if (PAR(cr))
814     {
815         /* After possible communicator splitting in make_dd_communicators.
816          * we can set up the intra/inter node communication.
817          */
818         gmx_setup_nodecomm(fplog, cr);
819     }
820
821     /* Initialize per-physical-node MPI process/thread ID and counters. */
822     gmx_init_intranode_counters(cr);
823 #if GMX_MPI
824     if (MULTISIM(cr))
825     {
826         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
827                 "This is simulation %d out of %d running as a composite GROMACS\n"
828                 "multi-simulation job. Setup for this simulation:\n",
829                 cr->ms->sim, cr->ms->nsim);
830     }
831     GMX_LOG(mdlog.warning).appendTextFormatted(
832             "Using %d MPI %s\n",
833             cr->nnodes,
834 #if GMX_THREAD_MPI
835             cr->nnodes == 1 ? "thread" : "threads"
836 #else
837             cr->nnodes == 1 ? "process" : "processes"
838 #endif
839             );
840     fflush(stderr);
841 #endif
842
843     /* Check and update hw_opt for the cut-off scheme */
844     check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
845
846     /* Check and update hw_opt for the number of MPI ranks */
847     check_and_update_hw_opt_3(&hw_opt);
848
849     gmx_omp_nthreads_init(mdlog, cr,
850                           hwinfo->nthreads_hw_avail,
851                           hw_opt.nthreads_omp,
852                           hw_opt.nthreads_omp_pme,
853                           (cr->duty & DUTY_PP) == 0,
854                           inputrec->cutoff_scheme == ecutsVERLET);
855
856 #ifndef NDEBUG
857     if (EI_TPI(inputrec->eI) &&
858         inputrec->cutoff_scheme == ecutsVERLET)
859     {
860         gmx_feenableexcept();
861     }
862 #endif
863
864     // Contains the ID of the GPU used by each PP rank on this node,
865     // indexed by that rank. Empty if no GPUs are selected for use on
866     // this node.
867     std::vector<int> gpuTaskAssignment;
868     if (nonbondedOnGpu)
869     {
870         /* Currently the DD code assigns duty to ranks that can
871          * include PP work that currently can be executed on a single
872          * GPU, if present and compatible.  This has to be coordinated
873          * across PP ranks on a node, with possible multiple devices
874          * or sharing devices on a node, either from the user
875          * selection, or automatically. */
876         bool rankCanUseGpu = cr->duty & DUTY_PP;
877         gpuTaskAssignment = mapPpRanksToGpus(rankCanUseGpu, cr, hwinfo->gpu_info, hw_opt);
878     }
879
880     reportGpuUsage(mdlog, hwinfo->gpu_info, !hw_opt.gpuIdTaskAssignment.empty(),
881                    gpuTaskAssignment, cr->nrank_pp_intranode, cr->nnodes > 1);
882
883     if (!gpuTaskAssignment.empty())
884     {
885         GMX_RELEASE_ASSERT(cr->nrank_pp_intranode == static_cast<int>(gpuTaskAssignment.size()),
886                            "The number of PP ranks on each node must equal the number of GPU tasks used on each node");
887     }
888
889     /* Prevent other ranks from continuing after an issue was found
890      * and reported as a fatal error.
891      *
892      * TODO This function implements a barrier so that MPI runtimes
893      * can organize an orderly shutdown if one of the ranks has had to
894      * issue a fatal error in various code already run. When we have
895      * MPI-aware error handling and reporting, this should be
896      * improved. */
897 #if GMX_MPI
898     if (PAR(cr))
899     {
900         MPI_Barrier(cr->mpi_comm_mysim);
901     }
902 #endif
903
904     /* Now that we know the setup is consistent, check for efficiency */
905     check_resource_division_efficiency(hwinfo, hw_opt.nthreads_tot, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
906                                        cr, mdlog);
907
908     gmx_device_info_t *shortRangedDeviceInfo = nullptr;
909     int                shortRangedDeviceId   = -1;
910     if (cr->duty & DUTY_PP)
911     {
912         if (!gpuTaskAssignment.empty())
913         {
914             shortRangedDeviceId   = gpuTaskAssignment[cr->rank_pp_intranode];
915             shortRangedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, shortRangedDeviceId);
916         }
917     }
918
919     if (DOMAINDECOMP(cr))
920     {
921         /* When we share GPUs over ranks, we need to know this for the DLB */
922         dd_setup_dlb_resource_sharing(cr, shortRangedDeviceId);
923     }
924
925     /* getting number of PP/PME threads
926        PME: env variable should be read only on one node to make sure it is
927        identical everywhere;
928      */
929     nthreads_pme = gmx_omp_nthreads_get(emntPME);
930
931     wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
932
933     if (PAR(cr))
934     {
935         /* Master synchronizes its value of reset_counters with all nodes
936          * including PME only nodes */
937         reset_counters = wcycle_get_reset_counters(wcycle);
938         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
939         wcycle_set_reset_counters(wcycle, reset_counters);
940     }
941
942     // Membrane embedding must be initialized before we call init_forcerec()
943     if (doMembed)
944     {
945         if (MASTER(cr))
946         {
947             fprintf(stderr, "Initializing membed");
948         }
949         /* Note that membed cannot work in parallel because mtop is
950          * changed here. Fix this if we ever want to make it run with
951          * multiple ranks. */
952         membed = init_membed(fplog, nfile, fnm, mtop, inputrec, globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
953     }
954
955     snew(nrnb, 1);
956     if (cr->duty & DUTY_PP)
957     {
958         /* Initiate forcerecord */
959         fr                 = mk_forcerec();
960         fr->forceProviders = mdModules.initForceProviders();
961         init_forcerec(fplog, mdlog, fr, fcd,
962                       inputrec, mtop, cr, box,
963                       opt2fn("-table", nfile, fnm),
964                       opt2fn("-tablep", nfile, fnm),
965                       getFilenm("-tableb", nfile, fnm),
966                       shortRangedDeviceInfo,
967                       FALSE,
968                       pforce);
969
970         /* Initialize QM-MM */
971         if (fr->bQMMM)
972         {
973             init_QMMMrec(cr, mtop, inputrec, fr);
974         }
975
976         /* Initialize the mdatoms structure.
977          * mdatoms is not filled with atom data,
978          * as this can not be done now with domain decomposition.
979          */
980         mdatoms = init_mdatoms(fplog, *mtop, *inputrec);
981
982         /* Initialize the virtual site communication */
983         vsite = init_vsite(mtop, cr, FALSE);
984
985         calc_shifts(box, fr->shift_vec);
986
987         /* With periodic molecules the charge groups should be whole at start up
988          * and the virtual sites should not be far from their proper positions.
989          */
990         if (!inputrec->bContinuation && MASTER(cr) &&
991             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
992         {
993             rvec *xGlobal = as_rvec_array(globalState->x.data());
994
995             /* Make molecules whole at start of run */
996             if (fr->ePBC != epbcNONE)
997             {
998                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, xGlobal);
999             }
1000             if (vsite)
1001             {
1002                 /* Correct initial vsite positions are required
1003                  * for the initial distribution in the domain decomposition
1004                  * and for the initial shell prediction.
1005                  */
1006                 construct_vsites_mtop(vsite, mtop, xGlobal);
1007             }
1008         }
1009
1010         if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1011         {
1012             ewaldcoeff_q  = fr->ic->ewaldcoeff_q;
1013             ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1014             pmedata       = &fr->pmedata;
1015         }
1016         else
1017         {
1018             pmedata = nullptr;
1019         }
1020     }
1021     else
1022     {
1023         /* This is a PME only node */
1024
1025         GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1026
1027         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1028         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1029         snew(pmedata, 1);
1030     }
1031
1032     if (hw_opt.thread_affinity != threadaffOFF)
1033     {
1034         /* Before setting affinity, check whether the affinity has changed
1035          * - which indicates that probably the OpenMP library has changed it
1036          * since we first checked).
1037          */
1038         gmx_check_thread_affinity_set(mdlog, cr,
1039                                       &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1040
1041         int nthread_local;
1042         /* threads on this MPI process or TMPI thread */
1043         if (cr->duty & DUTY_PP)
1044         {
1045             nthread_local = gmx_omp_nthreads_get(emntNonbonded);
1046         }
1047         else
1048         {
1049             nthread_local = gmx_omp_nthreads_get(emntPME);
1050         }
1051
1052         /* Set the CPU affinity */
1053         gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1054                                 nthread_local, nullptr);
1055     }
1056
1057     /* Initiate PME if necessary,
1058      * either on all nodes or on dedicated PME nodes only. */
1059     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1060     {
1061         if (mdatoms)
1062         {
1063             nChargePerturbed = mdatoms->nChargePerturbed;
1064             if (EVDW_PME(inputrec->vdwtype))
1065             {
1066                 nTypePerturbed   = mdatoms->nTypePerturbed;
1067             }
1068         }
1069         if (cr->npmenodes > 0)
1070         {
1071             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1072             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1073             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1074         }
1075
1076         if (cr->duty & DUTY_PME)
1077         {
1078             try
1079             {
1080                 gmx_device_info_t *pmeGpuInfo = nullptr;
1081                 auto               runMode    = PmeRunMode::CPU;
1082                 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1083                                       mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1084                                       mdrunOptions.reproducible,
1085                                       ewaldcoeff_q, ewaldcoeff_lj,
1086                                       nthreads_pme,
1087                                       runMode, nullptr, pmeGpuInfo, mdlog);
1088             }
1089             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1090             if (status != 0)
1091             {
1092                 gmx_fatal(FARGS, "Error %d initializing PME", status);
1093             }
1094         }
1095     }
1096
1097
1098     if (EI_DYNAMICS(inputrec->eI))
1099     {
1100         /* Turn on signal handling on all nodes */
1101         /*
1102          * (A user signal from the PME nodes (if any)
1103          * is communicated to the PP nodes.
1104          */
1105         signal_handler_install();
1106     }
1107
1108     if (cr->duty & DUTY_PP)
1109     {
1110         /* Assumes uniform use of the number of OpenMP threads */
1111         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1112
1113         if (inputrec->bPull)
1114         {
1115             /* Initialize pull code */
1116             inputrec->pull_work =
1117                 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1118                           mtop, cr, oenv, inputrec->fepvals->init_lambda,
1119                           EI_DYNAMICS(inputrec->eI) && MASTER(cr),
1120                           continuationOptions);
1121         }
1122
1123         if (inputrec->bRot)
1124         {
1125             /* Initialize enforced rotation code */
1126             init_rot(fplog, inputrec, nfile, fnm, cr, globalState.get(), mtop, oenv, mdrunOptions);
1127         }
1128
1129         /* Let init_constraints know whether we have essential dynamics constraints.
1130          * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1131          */
1132         bool doEdsam = (opt2fn_null("-ei", nfile, fnm) != nullptr || observablesHistory.edsamHistory);
1133
1134         constr = init_constraints(fplog, mtop, inputrec, doEdsam, cr);
1135
1136         if (DOMAINDECOMP(cr))
1137         {
1138             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1139             /* This call is not included in init_domain_decomposition mainly
1140              * because fr->cginfo_mb is set later.
1141              */
1142             dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1143                             domdecOptions.checkBondedInteractions,
1144                             fr->cginfo_mb);
1145         }
1146
1147         /* Now do whatever the user wants us to do (how flexible...) */
1148         my_integrator(inputrec->eI) (fplog, cr, mdlog, nfile, fnm,
1149                                      oenv,
1150                                      mdrunOptions,
1151                                      vsite, constr,
1152                                      mdModules.outputProvider(),
1153                                      inputrec, mtop,
1154                                      fcd,
1155                                      globalState.get(),
1156                                      &observablesHistory,
1157                                      mdatoms, nrnb, wcycle, fr,
1158                                      replExParams,
1159                                      membed,
1160                                      walltime_accounting);
1161
1162         if (inputrec->bRot)
1163         {
1164             finish_rot(inputrec->rot);
1165         }
1166
1167         if (inputrec->bPull)
1168         {
1169             finish_pull(inputrec->pull_work);
1170         }
1171
1172     }
1173     else
1174     {
1175         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1176         /* do PME only */
1177         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1178         gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec);
1179     }
1180
1181     wallcycle_stop(wcycle, ewcRUN);
1182
1183     /* Finish up, write some stuff
1184      * if rerunMD, don't write last frame again
1185      */
1186     finish_run(fplog, mdlog, cr,
1187                inputrec, nrnb, wcycle, walltime_accounting,
1188                fr ? fr->nbv : nullptr,
1189                fr ? fr->pmedata : nullptr,
1190                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1191
1192     // Free PME data
1193     if (pmedata)
1194     {
1195         gmx_pme_destroy(*pmedata); // TODO: pmedata is always a single element list, refactor
1196         pmedata = nullptr;
1197     }
1198
1199     /* Free GPU memory and context */
1200     free_gpu_resources(fr, cr, shortRangedDeviceInfo);
1201
1202     if (doMembed)
1203     {
1204         free_membed(membed);
1205     }
1206
1207     gmx_hardware_info_free(hwinfo);
1208
1209     /* Does what it says */
1210     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1211     walltime_accounting_destroy(walltime_accounting);
1212
1213     /* Close logfile already here if we were appending to it */
1214     if (MASTER(cr) && continuationOptions.appendFiles)
1215     {
1216         gmx_log_close(fplog);
1217     }
1218
1219     rc = (int)gmx_get_stop_condition();
1220
1221 #if GMX_THREAD_MPI
1222     /* we need to join all threads. The sub-threads join when they
1223        exit this function, but the master thread needs to be told to
1224        wait for that. */
1225     if (PAR(cr) && MASTER(cr))
1226     {
1227         tMPI_Finalize();
1228     }
1229 #endif
1230
1231     return rc;
1232 }
1233
1234 } // namespace gmx