2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * Implements functions in swapcoords.h.
39 * \author Carsten Kutzner <ckutzne@gwdg.de>
40 * \ingroup module_swap
44 #include "swapcoords.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/domdec/localatomset.h"
56 #include "gromacs/domdec/localatomsetmanager.h"
57 #include "gromacs/fileio/confio.h"
58 #include "gromacs/fileio/gmxfio.h"
59 #include "gromacs/fileio/xvgr.h"
60 #include "gromacs/gmxlib/network.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdlib/groupcoord.h"
63 #include "gromacs/mdrunutility/handlerestart.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/imdmodule.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/mdrunoptions.h"
69 #include "gromacs/mdtypes/observableshistory.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/mdtypes/swaphistory.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/timing/wallcycle.h"
74 #include "gromacs/topology/mtop_lookup.h"
75 #include "gromacs/topology/topology.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/pleasecite.h"
79 #include "gromacs/utility/smalloc.h"
80 #include "gromacs/utility/snprintf.h"
82 static const char* SwS = { "SWAP:" }; /**< For output that comes from the swap module */
83 static const char* SwSEmpty = { " " }; /**< Placeholder for multi-line output */
84 static const char* CompStr[eCompNR] = { "A", "B" }; /**< Compartment name */
85 static const char* SwapStr[eSwapTypesNR + 1] = { "", "X-", "Y-", "Z-",
86 nullptr }; /**< Name for the swap types. */
87 static const char* DimStr[DIM + 1] = { "X", "Y", "Z", nullptr }; /**< Name for the swap dimension. */
89 /** Keep track of through which channel the ions have passed */
97 static const char* ChannelString[eChHistNr] = { "none", "channel0", "channel1" }; /**< Name for the channels */
99 /*! \brief Domain identifier.
101 * Keeps track of from which compartment the ions came before passing the
111 static const char* DomainString[eDomainNr] = { "not_assigned", "Domain_A",
112 "Domain_B" }; /**< Name for the domains */
117 extern template LocalAtomSet LocalAtomSetManager::add<void, void>(ArrayRef<const int> globalAtomIndex);
120 * \brief Implement Computational Electrophysiology swapping.
122 class SwapCoordinates final : public IMDModule
125 IMdpOptionProvider* mdpOptionProvider() override { return nullptr; }
126 IMDOutputProvider* outputProvider() override { return nullptr; }
127 void initForceProviders(ForceProviders* /* forceProviders */) override {}
130 std::unique_ptr<IMDModule> createSwapCoordinatesModule()
132 return std::make_unique<SwapCoordinates>();
139 * Structure containing compartment-specific data.
141 typedef struct swap_compartment
143 int nMol; /**< Number of ion or water molecules detected
144 in this compartment. */
145 int nMolBefore; /**< Number of molecules before swapping. */
146 int nMolReq; /**< Requested number of molecules in compartment. */
147 real nMolAv; /**< Time-averaged number of molecules matching
148 the compartment conditions. */
149 int* nMolPast; /**< Past molecule counts for time-averaging. */
150 int* ind; /**< Indices to collective array of atoms. */
151 real* dist; /**< Distance of atom to bulk layer, which is
152 normally the center layer of the compartment */
153 int nalloc; /**< Allocation size for ind array. */
154 int inflow_net; /**< Net inflow of ions into this compartment. */
159 * This structure contains data needed for the groups involved in swapping:
160 * split group 0, split group 1, solvent group, ion groups.
162 typedef struct swap_group
164 /*!\brief Construct a swap group given the managed swap atoms.
166 * \param[in] atomset Managed indices of atoms that are part of the swap group.
168 swap_group(const gmx::LocalAtomSet& atomset);
169 char* molname = nullptr; /**< Name of the group or ion type */
170 int apm = 0; /**< Number of atoms in each molecule */
171 gmx::LocalAtomSet atomset; /**< The atom indices in the swap group */
172 rvec* xc = nullptr; /**< Collective array of group atom positions (size nat) */
173 ivec* xc_shifts = nullptr; /**< Current (collective) shifts (size nat) */
174 ivec* xc_eshifts = nullptr; /**< Extra shifts since last DD step (size nat) */
175 rvec* xc_old = nullptr; /**< Old (collective) positions (size nat) */
176 real q = 0.; /**< Total charge of one molecule of this group */
177 real* m = nullptr; /**< Masses (can be omitted, size apm) */
178 unsigned char* comp_from = nullptr; /**< (Collective) Stores from which compartment this
179 molecule has come. This way we keep track of
180 through which channel an ion permeates
181 (size nMol = nat/apm) */
182 unsigned char* comp_now = nullptr; /**< In which compartment this ion is now (size nMol) */
183 unsigned char* channel_label = nullptr; /**< Which channel was passed at last by this ion?
185 rvec center; /**< Center of the group; COM if masses are used */
186 t_compartment comp[eCompNR]; /**< Distribution of particles of this group across
187 the two compartments */
188 real vacancy[eCompNR]; /**< How many molecules need to be swapped in? */
189 int fluxfromAtoB[eChanNR]; /**< Net flux of ions per channel */
190 int nCyl[eChanNR]; /**< Number of ions residing in a channel */
191 int nCylBoth = 0; /**< Ions assigned to cyl0 and cyl1. Not good. */
194 t_swapgrp::swap_group(const gmx::LocalAtomSet& atomset) : atomset{ atomset }
199 for (int compartment = eCompA; compartment < eCompNR; ++compartment)
201 comp[compartment] = {};
202 vacancy[compartment] = 0;
204 for (int channel = eChan0; channel < eChanNR; ++channel)
206 fluxfromAtoB[channel] = 0;
212 * Main (private) data structure for the position swapping protocol.
216 int swapdim; /**< One of XX, YY, ZZ */
217 t_pbc* pbc; /**< Needed to make molecules whole. */
218 FILE* fpout; /**< Output file. */
219 int ngrp; /**< Number of t_swapgrp groups */
220 std::vector<t_swapgrp> group; /**< Separate groups for channels, solvent, ions */
221 int fluxleak; /**< Flux not going through any of the channels. */
222 real deltaQ; /**< The charge imbalance between the compartments. */
226 /*! \brief Check whether point is in channel.
228 * A channel is a cylinder defined by a disc
229 * with radius r around its center c. The thickness of the cylinder is
236 * <---------c--------->
242 * \param[in] point The position (xyz) under consideration.
243 * \param[in] center The center of the cylinder.
244 * \param[in] d_up The upper extension of the cylinder.
245 * \param[in] d_down The lower extension.
246 * \param[in] r_cyl2 Cylinder radius squared.
247 * \param[in] pbc Structure with info about periodic boundary conditions.
248 * \param[in] normal The membrane normal direction is typically 3, i.e. z, but can be x or y also.
250 * \returns Whether the point is inside the defined cylindric channel.
252 static gmx_bool is_in_channel(rvec point, rvec center, real d_up, real d_down, real r_cyl2, t_pbc* pbc, int normal)
255 int plane1, plane2; /* Directions tangential to membrane */
258 plane1 = (normal + 1) % 3; /* typically 0, i.e. XX */
259 plane2 = (normal + 2) % 3; /* typically 1, i.e. YY */
261 /* Get the distance vector dr between the point and the center of the cylinder */
262 pbc_dx(pbc, point, center, dr); /* This puts center in the origin */
264 /* Check vertical direction */
265 if ((dr[normal] > d_up) || (dr[normal] < -d_down))
270 /* Check radial direction */
271 if ((dr[plane1] * dr[plane1] + dr[plane2] * dr[plane2]) > r_cyl2)
276 /* All check passed, this point is in the cylinder */
281 /*! \brief Prints output to CompEL output file.
283 * Prints to swap output file how many ions are in each compartment,
284 * where the centers of the split groups are, and how many ions of each type
285 * passed the channels.
287 static void print_ionlist(t_swap* s, double time, const char comment[])
290 fprintf(s->fpout, "%12.5e", time);
292 // Output number of molecules and difference to reference counts for each
293 // compartment and ion type
294 for (int iComp = 0; iComp < eCompNR; iComp++)
296 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
298 t_compartment* comp = &s->group[ig].comp[iComp];
300 fprintf(s->fpout, "%10d%10.1f%10d", comp->nMol, comp->nMolAv - comp->nMolReq, comp->inflow_net);
304 // Output center of split groups
305 fprintf(s->fpout, "%10g%10g", s->group[eGrpSplit0].center[s->swapdim],
306 s->group[eGrpSplit1].center[s->swapdim]);
308 // Output ion flux for each channel and ion type
309 for (int iChan = 0; iChan < eChanNR; iChan++)
311 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
313 t_swapgrp* g = &s->group[ig];
314 fprintf(s->fpout, "%10d", g->fluxfromAtoB[iChan]);
318 /* Output the number of molecules that leaked from A to B */
319 fprintf(s->fpout, "%10d", s->fluxleak);
321 fprintf(s->fpout, "%s\n", comment);
325 /*! \brief Get the center of a group of nat atoms.
327 * Since with PBC an atom group might not be whole, use the first atom as the
328 * reference atom and determine the center with respect to this reference.
330 static void get_molecule_center(rvec x[], int nat, const real* weights, rvec center, t_pbc* pbc)
333 rvec weightedPBCimage;
335 rvec reference, correctPBCimage, dx;
338 /* Use the first atom as the reference and put other atoms near that one */
339 /* This does not work for large molecules that span > half of the box! */
340 copy_rvec(x[0], reference);
342 /* Calculate either the weighted center or simply the center of geometry */
345 for (i = 0; i < nat; i++)
347 /* PBC distance between position and reference */
348 pbc_dx(pbc, x[i], reference, dx);
350 /* Add PBC distance to reference */
351 rvec_add(reference, dx, correctPBCimage);
353 /* Take weight into account */
354 if (nullptr == weights)
363 svmul(wi, correctPBCimage, weightedPBCimage);
366 rvec_inc(center, weightedPBCimage);
370 svmul(1.0 / wsum, center, center);
374 /*! \brief Return TRUE if position x of ion (or water) is found in the compartment,
375 * i.e. between w1 and w2.
377 * One can define and additional offset "b" if one wants to exchange ions/water
378 * to or from a plane not directly in the middle of w1 and w2. The offset can be
379 * in ]-1.0, ..., +1.0 [.
380 * A bulkOffset of 0.0 means 'no offset', so the swap-layer is directly in the
381 * middle between w1 and w2. Offsets -1.0 < b < 0.0 will yield swaps nearer to w1,
382 * whereas offsets 0.0 < 0 < +1.0 will yield swaps nearer to w2.
386 * ||--------------+-------------|-------------+------------------------||
387 * w1 ? ? ? ? ? ? ? ? ? ? ? w2
388 * ||--------------+-------------|----b--------+------------------------||
393 * \param[in] w1 Position of 'wall' atom 1.
394 * \param[in] w2 Position of 'wall' atom 2.
395 * \param[in] x Position of the ion or the water molecule under consideration.
396 * \param[in] l Length of the box, from || to || in the sketch.
397 * \param[in] bulkOffset Where is the bulk layer "b" to be found between w1 and w2?
398 * \param[out] distance_from_b Distance of x to the bulk layer "b".
400 * \returns TRUE if x is between w1 and w2.
402 * Also computes the distance of x to the compartment center (the layer that is
403 * normally situated in the middle of w1 and w2 that would be considered as having
404 * the bulk concentration of ions).
406 static gmx_bool compartment_contains_atom(real w1, real w2, real x, real l, real bulkOffset, real* distance_from_b)
412 /* First set the origin in the middle of w1 and w2 */
419 /* Now choose the PBC image of x that is closest to the origin: */
430 *distance_from_b = static_cast<real>(fabs(x - bulkOffset * 0.5 * width));
432 /* Return TRUE if we now are in area "????" */
433 return (x >= w1) && (x < w2);
437 /*! \brief Updates the time-averaged number of ions in a compartment. */
438 static void update_time_window(t_compartment* comp, int values, int replace)
444 /* Put in the new value */
447 comp->nMolPast[replace] = comp->nMol;
450 /* Compute the new time-average */
452 for (i = 0; i < values; i++)
454 average += comp->nMolPast[i];
457 comp->nMolAv = average;
461 /*! \brief Add the atom with collective index ci to the atom list in compartment 'comp'.
463 * \param[in] ci Index of this ion in the collective xc array.
464 * \param[inout] comp Compartment to add this atom to.
465 * \param[in] distance Shortest distance of this atom to the bulk layer,
466 * from which ion/water pairs are selected for swapping.
468 static void add_to_list(int ci, t_compartment* comp, real distance)
472 if (nr >= comp->nalloc)
474 comp->nalloc = over_alloc_dd(nr + 1);
475 srenew(comp->ind, comp->nalloc);
476 srenew(comp->dist, comp->nalloc);
479 comp->dist[nr] = distance;
484 /*! \brief Determine the compartment boundaries from the channel centers. */
485 static void get_compartment_boundaries(int c, t_swap* s, const matrix box, real* left, real* right)
488 real leftpos, rightpos, leftpos_orig;
493 gmx_fatal(FARGS, "No compartment %c.", c + 'A');
496 pos0 = s->group[eGrpSplit0].center[s->swapdim];
497 pos1 = s->group[eGrpSplit1].center[s->swapdim];
510 /* This gets us the other compartment: */
513 leftpos_orig = leftpos;
515 rightpos = leftpos_orig + box[s->swapdim][s->swapdim];
523 /*! \brief Determine the per-channel ion flux.
525 * To determine the flux through the individual channels, we
526 * remember the compartment and channel history of each ion. An ion can be
527 * either in channel0 or channel1, or in the remaining volume of compartment
531 * +-----------------+
534 * ||||||||||0|||||||| bilayer with channel 0
539 * |||||1||||||||||||| bilayer with channel 1
542 * +-----------------+
546 static void detect_flux_per_channel(t_swapgrp* g,
550 unsigned char* comp_now,
551 unsigned char* comp_from,
552 unsigned char* channel_label,
562 gmx_bool in_cyl0, in_cyl1;
568 /* Check whether ion is inside any of the channels */
569 in_cyl0 = is_in_channel(atomPosition, s->group[eGrpSplit0].center, sc->cyl0u, sc->cyl0l,
570 cyl0_r2, s->pbc, sd);
571 in_cyl1 = is_in_channel(atomPosition, s->group[eGrpSplit1].center, sc->cyl1u, sc->cyl1l,
572 cyl1_r2, s->pbc, sd);
574 if (in_cyl0 && in_cyl1)
576 /* Ion appears to be in both channels. Something is severely wrong! */
578 *comp_now = eDomainNotset;
579 *comp_from = eDomainNotset;
580 *channel_label = eChHistPassedNone;
584 /* Ion is in channel 0 now */
585 *channel_label = eChHistPassedCh0;
586 *comp_now = eDomainNotset;
591 /* Ion is in channel 1 now */
592 *channel_label = eChHistPassedCh1;
593 *comp_now = eDomainNotset;
598 /* Ion is not in any of the channels, so it must be in domain A or B */
601 *comp_now = eDomainA;
605 *comp_now = eDomainB;
609 /* Only take action, if ion is now in domain A or B, and was before
610 * in the other domain!
612 if (eDomainNotset == *comp_from)
614 /* Maybe we can set the domain now */
615 *comp_from = *comp_now; /* Could still be eDomainNotset, though */
617 else if ((*comp_now != eDomainNotset) /* if in channel */
618 && (*comp_from != *comp_now))
620 /* Obviously the ion changed its domain.
621 * Count this for the channel through which it has passed. */
622 switch (*channel_label)
624 case eChHistPassedNone:
627 fprintf(stderr, " %s Warning! Step %s, ion %d moved from %s to %s\n", SwS,
628 gmx_step_str(step, buf), iAtom, DomainString[*comp_from],
629 DomainString[*comp_now]);
632 fprintf(stderr, ", possibly due to a swap in the original simulation.\n");
637 "but did not pass cyl0 or cyl1 as defined in the .mdp file.\n"
638 "Do you have an ion somewhere within the membrane?\n");
639 /* Write this info to the CompEL output file: */
641 " # Warning: step %s, ion %d moved from %s to %s (probably through the "
643 gmx_step_str(step, buf), iAtom, DomainString[*comp_from],
644 DomainString[*comp_now]);
647 case eChHistPassedCh0:
648 case eChHistPassedCh1:
649 if (*channel_label == eChHistPassedCh0)
658 if (eDomainA == *comp_from)
660 g->fluxfromAtoB[chan_nr]++;
664 g->fluxfromAtoB[chan_nr]--;
666 fprintf(fpout, "# Atom nr. %d finished passing %s.\n", iAtom,
667 ChannelString[*channel_label]);
670 gmx_fatal(FARGS, "%s Unknown channel history entry for ion type '%s'\n", SwS, g->molname);
673 /* This ion has moved to the _other_ compartment ... */
674 *comp_from = *comp_now;
675 /* ... and it did not pass any channel yet */
676 *channel_label = eChHistPassedNone;
681 /*! \brief Determines which ions or solvent molecules are in compartment A and B */
682 static void sortMoleculesIntoCompartments(t_swapgrp* g,
692 int nMolNotInComp[eCompNR]; /* consistency check */
693 real cyl0_r2 = sc->cyl0r * sc->cyl0r;
694 real cyl1_r2 = sc->cyl1r * sc->cyl1r;
696 /* Get us a counter that cycles in the range of [0 ... sc->nAverage[ */
697 int replace = (step / sc->nstswap) % sc->nAverage;
699 for (int comp = eCompA; comp <= eCompB; comp++)
703 /* Get lists of atoms that match criteria for this compartment */
704 get_compartment_boundaries(comp, s, box, &left, &right);
706 /* First clear the ion molecule lists */
707 g->comp[comp].nMol = 0;
708 nMolNotInComp[comp] = 0; /* consistency check */
710 /* Loop over the molecules and atoms of this group */
711 for (int iMol = 0, iAtom = 0; iAtom < static_cast<int>(g->atomset.numAtomsGlobal());
712 iAtom += g->apm, iMol++)
717 /* Is this first atom of the molecule in the compartment that we look at? */
718 if (compartment_contains_atom(left, right, g->xc[iAtom][sd], box[sd][sd],
719 sc->bulkOffset[comp], &dist))
721 /* Add the first atom of this molecule to the list of molecules in this compartment */
722 add_to_list(iAtom, &g->comp[comp], dist);
724 /* Master also checks for ion groups through which channel each ion has passed */
725 if (MASTER(cr) && (g->comp_now != nullptr) && !bIsSolvent)
727 int globalAtomNr = g->atomset.globalIndex()[iAtom] + 1; /* PDB index starts at 1 ... */
728 detect_flux_per_channel(g, globalAtomNr, comp, g->xc[iAtom], &g->comp_now[iMol],
729 &g->comp_from[iMol], &g->channel_label[iMol], sc, s,
730 cyl0_r2, cyl1_r2, step, bRerun, fpout);
735 nMolNotInComp[comp]++;
738 /* Correct the time-averaged number of ions in the compartment */
741 update_time_window(&g->comp[comp], sc->nAverage, replace);
745 /* Flux detection warnings */
746 if (MASTER(cr) && !bIsSolvent)
752 "%s Warning: %d atoms were detected as being in both channels! Probably your "
754 "%s cylinder is way too large, or one compartment has collapsed (step "
756 SwS, g->nCylBoth, SwS, step);
758 fprintf(s->fpout, "Warning: %d atoms were assigned to both channels!\n", g->nCylBoth);
764 if (bIsSolvent && nullptr != fpout)
766 fprintf(fpout, "# Solv. molecules in comp.%s: %d comp.%s: %d\n", CompStr[eCompA],
767 g->comp[eCompA].nMol, CompStr[eCompB], g->comp[eCompB].nMol);
770 /* Consistency checks */
771 const auto numMolecules = static_cast<int>(g->atomset.numAtomsGlobal() / g->apm);
772 if (nMolNotInComp[eCompA] + nMolNotInComp[eCompB] != numMolecules)
775 "%s Warning: Inconsistency while assigning '%s' molecules to compartments. !inA: "
776 "%d, !inB: %d, total molecules %d\n",
777 SwS, g->molname, nMolNotInComp[eCompA], nMolNotInComp[eCompB], numMolecules);
780 int sum = g->comp[eCompA].nMol + g->comp[eCompB].nMol;
781 if (sum != numMolecules)
784 "%s Warning: %d molecules are in group '%s', but altogether %d have been assigned "
785 "to the compartments.\n",
786 SwS, numMolecules, g->molname, sum);
791 /*! \brief Find out how many group atoms are in the compartments initially */
792 static void get_initial_ioncounts(const t_inputrec* ir,
794 const rvec x[], /* the initial positions */
804 /* Loop over the user-defined (ion) groups */
805 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
809 /* Copy the initial positions of the atoms in the group
810 * to the collective array so that we can compartmentalize */
811 for (size_t i = 0; i < g->atomset.numAtomsGlobal(); i++)
813 int ind = g->atomset.globalIndex()[i];
814 copy_rvec(x[ind], g->xc[i]);
817 /* Set up the compartments and get lists of atoms in each compartment */
818 sortMoleculesIntoCompartments(g, cr, sc, s, box, 0, s->fpout, bRerun, FALSE);
820 /* Set initial molecule counts if requested (as signaled by "-1" value) */
821 for (int ic = 0; ic < eCompNR; ic++)
823 int requested = sc->grp[ig].nmolReq[ic];
826 g->comp[ic].nMolReq = g->comp[ic].nMol;
830 g->comp[ic].nMolReq = requested;
834 /* Check whether the number of requested molecules adds up to the total number */
835 int req = g->comp[eCompA].nMolReq + g->comp[eCompB].nMolReq;
836 int tot = g->comp[eCompA].nMol + g->comp[eCompB].nMol;
841 "Mismatch of the number of %s ions summed over both compartments.\n"
842 "You requested a total of %d ions (%d in A and %d in B),\n"
843 "but there are a total of %d ions of this type in the system.\n",
844 g->molname, req, g->comp[eCompA].nMolReq, g->comp[eCompB].nMolReq, tot);
847 /* Initialize time-averaging:
848 * Write initial concentrations to all time bins to start with */
849 for (int ic = 0; ic < eCompNR; ic++)
851 g->comp[ic].nMolAv = g->comp[ic].nMol;
852 for (int i = 0; i < sc->nAverage; i++)
854 g->comp[ic].nMolPast[i] = g->comp[ic].nMol;
861 /*! \brief Copy history of ion counts from checkpoint file.
863 * When called, the checkpoint file has already been read in. Here we copy
864 * over the values from .cpt file to the swap data structure.
866 static void get_initial_ioncounts_from_cpt(const t_inputrec* ir,
868 swaphistory_t* swapstate,
880 /* Copy the past values from the checkpoint values that have been read in already */
883 fprintf(stderr, "%s Copying values from checkpoint\n", SwS);
886 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
889 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
891 for (int ic = 0; ic < eCompNR; ic++)
893 g->comp[ic].nMolReq = gs->nMolReq[ic];
894 g->comp[ic].inflow_net = gs->inflow_net[ic];
898 fprintf(stderr, "%s ... Influx netto: %d Requested: %d Past values: ", SwS,
899 g->comp[ic].inflow_net, g->comp[ic].nMolReq);
902 for (int j = 0; j < sc->nAverage; j++)
904 g->comp[ic].nMolPast[j] = gs->nMolPast[ic][j];
907 fprintf(stderr, "%d ", g->comp[ic].nMolPast[j]);
912 fprintf(stderr, "\n");
920 /*! \brief The master lets all others know about the initial ion counts. */
921 static void bc_initial_concentrations(t_commrec* cr, t_swapcoords* swap, t_swap* s)
927 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
931 for (ic = 0; ic < eCompNR; ic++)
933 gmx_bcast(sizeof(g->comp[ic].nMolReq), &(g->comp[ic].nMolReq), cr);
934 gmx_bcast(sizeof(g->comp[ic].nMol), &(g->comp[ic].nMol), cr);
935 gmx_bcast(swap->nAverage * sizeof(g->comp[ic].nMolPast[0]), g->comp[ic].nMolPast, cr);
941 /*! \brief Ensure that each atom belongs to at most one of the swap groups. */
942 static void check_swap_groups(t_swap* s, int nat, gmx_bool bVerbose)
944 int* nGroup = nullptr; /* This array counts for each atom in the MD system to
945 how many swap groups it belongs (should be 0 or 1!) */
947 int nMultiple = 0; /* Number of atoms belonging to multiple groups */
952 fprintf(stderr, "%s Making sure each atom belongs to at most one of the swap groups.\n", SwS);
955 /* Add one to the group count of atoms belonging to a swap group: */
957 for (int i = 0; i < s->ngrp; i++)
959 t_swapgrp* g = &s->group[i];
960 for (size_t j = 0; j < g->atomset.numAtomsGlobal(); j++)
962 /* Get the global index of this atom of this group: */
963 ind = g->atomset.globalIndex()[j];
967 /* Make sure each atom belongs to at most one of the groups: */
968 for (int i = 0; i < nat; i++)
980 "%s Cannot perform swapping since %d atom%s allocated to more than one swap "
982 "%s Each atom must be allocated to at most one of the split groups, the swap "
983 "groups, or the solvent.\n"
984 "%s Check the .mdp file settings regarding the swap index groups or the index "
985 "groups themselves.\n",
986 SwS, nMultiple, (1 == nMultiple) ? " is" : "s are", SwSEmpty, SwSEmpty);
991 /*! \brief Get the number of atoms per molecule for this group.
993 * Also ensure that all the molecules in this group have this number of atoms.
995 static int get_group_apm_check(int igroup, t_swap* s, gmx_bool bVerbose, gmx_mtop_t* mtop)
997 t_swapgrp* g = &s->group[igroup];
998 const int* ind = s->group[igroup].atomset.globalIndex().data();
999 int nat = s->group[igroup].atomset.numAtomsGlobal();
1001 /* Determine the number of solvent atoms per solvent molecule from the
1002 * first solvent atom: */
1004 mtopGetMolblockIndex(mtop, ind[0], &molb, nullptr, nullptr);
1005 int apm = mtop->moleculeBlockIndices[molb].numAtomsPerMolecule;
1009 fprintf(stderr, "%s Checking whether all %s molecules consist of %d atom%s\n", SwS,
1010 g->molname, apm, apm > 1 ? "s" : "");
1013 /* Check whether this is also true for all other solvent atoms */
1014 for (int i = 1; i < nat; i++)
1016 mtopGetMolblockIndex(mtop, ind[i], &molb, nullptr, nullptr);
1017 if (apm != mtop->moleculeBlockIndices[molb].numAtomsPerMolecule)
1019 gmx_fatal(FARGS, "Not all molecules of swap group %d consist of %d atoms.", igroup, apm);
1023 // TODO: check whether charges and masses of each molecule are identical!
1028 /*! \brief Print the legend to the swap output file.
1030 * Also print the initial values of ion counts and position of split groups.
1032 static void print_ionlist_legend(const t_inputrec* ir, t_swap* s, const gmx_output_env_t* oenv)
1034 const char** legend;
1038 int nIonTypes = ir->swap->ngrp - eSwapFixedGrpNR;
1039 snew(legend, eCompNR * nIonTypes * 3 + 2 + eChanNR * nIonTypes + 1);
1041 // Number of molecules and difference to reference counts for each
1042 // compartment and ion type
1043 for (int ic = count = 0; ic < eCompNR; ic++)
1045 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1047 t_swapGroup* g = &ir->swap->grp[ig];
1048 real q = s->group[ig].q;
1050 snprintf(buf, STRLEN, "%s %s ions (charge %s%g)", CompStr[ic], g->molname, q > 0 ? "+" : "", q);
1051 legend[count++] = gmx_strdup(buf);
1053 snprintf(buf, STRLEN, "%s av. mismatch to %d %s ions", CompStr[ic],
1054 s->group[ig].comp[ic].nMolReq, g->molname);
1055 legend[count++] = gmx_strdup(buf);
1057 snprintf(buf, STRLEN, "%s net %s ion influx", CompStr[ic], g->molname);
1058 legend[count++] = gmx_strdup(buf);
1062 // Center of split groups
1063 snprintf(buf, STRLEN, "%scenter of %s of split group 0", SwapStr[ir->eSwapCoords],
1064 (nullptr != s->group[eGrpSplit0].m) ? "mass" : "geometry");
1065 legend[count++] = gmx_strdup(buf);
1066 snprintf(buf, STRLEN, "%scenter of %s of split group 1", SwapStr[ir->eSwapCoords],
1067 (nullptr != s->group[eGrpSplit1].m) ? "mass" : "geometry");
1068 legend[count++] = gmx_strdup(buf);
1070 // Ion flux for each channel and ion type
1071 for (int ic = 0; ic < eChanNR; ic++)
1073 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1075 t_swapGroup* g = &ir->swap->grp[ig];
1076 snprintf(buf, STRLEN, "A->ch%d->B %s permeations", ic, g->molname);
1077 legend[count++] = gmx_strdup(buf);
1081 // Number of molecules that leaked from A to B
1082 snprintf(buf, STRLEN, "leakage");
1083 legend[count++] = gmx_strdup(buf);
1085 xvgr_legend(s->fpout, count, legend, oenv);
1088 "# Instantaneous ion counts and time-averaged differences to requested numbers\n");
1090 // We add a simple text legend helping to identify the columns with xvgr legend strings
1091 fprintf(s->fpout, "# time (ps)");
1092 for (int i = 0; i < count; i++)
1094 snprintf(buf, STRLEN, "s%d", i);
1095 fprintf(s->fpout, "%10s", buf);
1097 fprintf(s->fpout, "\n");
1102 /*! \brief Initialize channel ion flux detection routine.
1104 * Initialize arrays that keep track of where the ions come from and where
1107 static void detect_flux_per_channel_init(t_swap* s, swaphistory_t* swapstate, const bool isRestart)
1110 swapstateIons_t* gs;
1112 /* All these flux detection routines run on the master only */
1113 if (swapstate == nullptr)
1118 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1121 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1123 /******************************************************/
1124 /* Channel and domain history for the individual ions */
1125 /******************************************************/
1126 if (isRestart) /* set the pointers right */
1128 g->comp_from = gs->comp_from;
1129 g->channel_label = gs->channel_label;
1131 else /* allocate memory for molecule counts */
1133 snew(g->comp_from, g->atomset.numAtomsGlobal() / g->apm);
1134 gs->comp_from = g->comp_from;
1135 snew(g->channel_label, g->atomset.numAtomsGlobal() / g->apm);
1136 gs->channel_label = g->channel_label;
1138 snew(g->comp_now, g->atomset.numAtomsGlobal() / g->apm);
1140 /* Initialize the channel and domain history counters */
1141 for (size_t i = 0; i < g->atomset.numAtomsGlobal() / g->apm; i++)
1143 g->comp_now[i] = eDomainNotset;
1146 g->comp_from[i] = eDomainNotset;
1147 g->channel_label[i] = eChHistPassedNone;
1151 /************************************/
1152 /* Channel fluxes for both channels */
1153 /************************************/
1154 g->nCyl[eChan0] = 0;
1155 g->nCyl[eChan1] = 0;
1161 fprintf(stderr, "%s Copying channel fluxes from checkpoint file data\n", SwS);
1165 // Loop over ion types (and both channels)
1166 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1169 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1171 for (int ic = 0; ic < eChanNR; ic++)
1173 fprintf(stderr, "%s Channel %d flux history for ion type %s (charge %g): ", SwS, ic,
1177 g->fluxfromAtoB[ic] = gs->fluxfromAtoB[ic];
1181 g->fluxfromAtoB[ic] = 0;
1184 fprintf(stderr, "%d molecule%s", g->fluxfromAtoB[ic], g->fluxfromAtoB[ic] == 1 ? "" : "s");
1185 fprintf(stderr, "\n");
1189 /* Set pointers for checkpoint writing */
1190 swapstate->fluxleak_p = &s->fluxleak;
1191 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1194 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1196 for (int ic = 0; ic < eChanNR; ic++)
1198 gs->fluxfromAtoB_p[ic] = &g->fluxfromAtoB[ic];
1204 /*! \brief Outputs the initial structure to PDB file for debugging reasons.
1206 * Output the starting structure so that in case of multimeric channels
1207 * the user can check whether we have the correct PBC image for all atoms.
1208 * If this is not correct, the ion counts per channel will be very likely
1211 static void outputStartStructureIfWanted(gmx_mtop_t* mtop, rvec* x, int ePBC, const matrix box)
1213 char* env = getenv("GMX_COMPELDUMP");
1218 "\n%s Found env.var. GMX_COMPELDUMP, will output CompEL starting structure made "
1220 "%s In case of multimeric channels, please check whether they have the correct PBC "
1221 "representation.\n",
1224 write_sto_conf_mtop("CompELAssumedWholeConfiguration.pdb", *mtop->name, mtop, x, nullptr,
1230 /*! \brief Initialize the swapstate structure, used for checkpoint writing.
1232 * The swapstate struct stores the information we need to make the channels
1233 * whole again after restarts from a checkpoint file. Here we do the following:
1234 * a) If we did not start from .cpt, we prepare the struct for proper .cpt writing,
1235 * b) if we did start from .cpt, we copy over the last whole structures from .cpt,
1236 * c) in any case, for subsequent checkpoint writing, we set the pointers in
1237 * swapstate to the x_old arrays, which contain the correct PBC representation of
1238 * multimeric channels at the last time step.
1240 static void init_swapstate(swaphistory_t* swapstate,
1244 const rvec* x, /* the initial positions */
1246 const t_inputrec* ir)
1248 rvec* x_pbc = nullptr; /* positions of the whole MD system with molecules made whole */
1252 /* We always need the last whole positions such that
1253 * in the next time step we can make the channels whole again in PBC */
1254 if (swapstate->bFromCpt)
1256 /* Copy the last whole positions of each channel from .cpt */
1257 g = &(s->group[eGrpSplit0]);
1258 for (size_t i = 0; i < g->atomset.numAtomsGlobal(); i++)
1260 copy_rvec(swapstate->xc_old_whole[eChan0][i], g->xc_old[i]);
1262 g = &(s->group[eGrpSplit1]);
1263 for (size_t i = 0; i < g->atomset.numAtomsGlobal(); i++)
1265 copy_rvec(swapstate->xc_old_whole[eChan1][i], g->xc_old[i]);
1270 swapstate->eSwapCoords = ir->eSwapCoords;
1272 /* Set the number of ion types and allocate memory for checkpointing */
1273 swapstate->nIonTypes = s->ngrp - eSwapFixedGrpNR;
1274 snew(swapstate->ionType, swapstate->nIonTypes);
1276 /* Store the total number of ions of each type in the swapstateIons
1277 * structure that is accessible during checkpoint writing */
1278 for (int ii = 0; ii < swapstate->nIonTypes; ii++)
1280 swapstateIons_t* gs = &swapstate->ionType[ii];
1281 gs->nMol = sc->grp[ii + eSwapFixedGrpNR].nat;
1284 /* Extract the initial split group positions. */
1286 /* Remove pbc, make molecule whole. */
1287 snew(x_pbc, mtop->natoms);
1288 copy_rvecn(x, x_pbc, 0, mtop->natoms);
1290 /* This can only make individual molecules whole, not multimers */
1291 do_pbc_mtop(ir->ePBC, box, mtop, x_pbc);
1293 /* Output the starting structure? */
1294 outputStartStructureIfWanted(mtop, x_pbc, ir->ePBC, box);
1296 /* If this is the first run (i.e. no checkpoint present) we assume
1297 * that the starting positions give us the correct PBC representation */
1298 for (int ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1300 g = &(s->group[ig]);
1301 for (size_t i = 0; i < g->atomset.numAtomsGlobal(); i++)
1303 copy_rvec(x_pbc[g->atomset.globalIndex()[i]], g->xc_old[i]);
1308 /* Prepare swapstate arrays for later checkpoint writing */
1309 swapstate->nat[eChan0] = s->group[eGrpSplit0].atomset.numAtomsGlobal();
1310 swapstate->nat[eChan1] = s->group[eGrpSplit1].atomset.numAtomsGlobal();
1313 /* For subsequent checkpoint writing, set the swapstate pointers to the xc_old
1314 * arrays that get updated at every swapping step */
1315 swapstate->xc_old_whole_p[eChan0] = &s->group[eGrpSplit0].xc_old;
1316 swapstate->xc_old_whole_p[eChan1] = &s->group[eGrpSplit1].xc_old;
1319 /*! \brief Determine the total charge imbalance resulting from the swap groups */
1320 static real getRequestedChargeImbalance(t_swap* s)
1325 real particle_charge;
1326 real particle_number[eCompNR];
1328 // s->deltaQ = ( (-1) * s->comp[eCompA][eIonNEG].nat_req + s->comp[eCompA][eIonPOS].nat_req )
1329 // - ( (-1) * s->comp[eCompB][eIonNEG].nat_req + s->comp[eCompB][eIonPOS].nat_req );
1331 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1335 particle_charge = g->q;
1336 particle_number[eCompA] = g->comp[eCompA].nMolReq;
1337 particle_number[eCompB] = g->comp[eCompB].nMolReq;
1339 DeltaQ += particle_charge * (particle_number[eCompA] - particle_number[eCompB]);
1346 /*! \brief Sorts anions and cations into two separate groups
1348 * This routine should be called for the 'anions' and 'cations' group,
1349 * of which the indices were lumped together in the older version of the code.
1351 static void copyIndicesToGroup(const int* indIons, int nIons, t_swapGroup* g, t_commrec* cr)
1355 /* If explicit ion counts were requested in the .mdp file
1356 * (by setting positive values for the number of ions),
1357 * we can make an additional consistency check here */
1358 if ((g->nmolReq[eCompA] < 0) && (g->nmolReq[eCompB] < 0))
1360 if (g->nat != (g->nmolReq[eCompA] + g->nmolReq[eCompB]))
1362 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
1363 "%s Inconsistency while importing swap-related data from an old "
1364 "input file version.\n"
1365 "%s The requested ion counts in compartments A (%d) and B (%d)\n"
1366 "%s do not add up to the number of ions (%d) of this type for the "
1368 SwS, SwSEmpty, g->nmolReq[eCompA], g->nmolReq[eCompB], SwSEmpty,
1369 g->nat, g->molname);
1373 srenew(g->ind, g->nat);
1374 for (int i = 0; i < g->nat; i++)
1376 g->ind[i] = indIons[i];
1381 /*! \brief Converts old .tpr file CompEL contents to new data layout.
1383 * If we have read an old .tpr file (tpxv <= tpxv_CompElPolyatomicIonsAndMultipleIonTypes),
1384 * anions and cations are stored together in group #3. In the new
1385 * format we store each ion type in a separate group.
1386 * The 'classic' groups are:
1387 * #0 split group 0 - OK
1388 * #1 split group 1 - OK
1390 * #3 anions - contains also cations, needs to be converted
1391 * #4 cations - empty before conversion
1394 static void convertOldToNewGroupFormat(t_swapcoords* sc, gmx_mtop_t* mtop, gmx_bool bVerbose, t_commrec* cr)
1396 t_swapGroup* g = &sc->grp[3];
1398 /* Loop through the atom indices of group #3 (anions) and put all indices
1399 * that belong to cations into the cation group.
1403 int* indAnions = nullptr;
1404 int* indCations = nullptr;
1405 snew(indAnions, g->nat);
1406 snew(indCations, g->nat);
1409 for (int i = 0; i < g->nat; i++)
1411 const t_atom& atom = mtopGetAtomParameters(mtop, g->ind[i], &molb);
1414 // This is an anion, add it to the list of anions
1415 indAnions[nAnions++] = g->ind[i];
1419 // This is a cation, add it to the list of cations
1420 indCations[nCations++] = g->ind[i];
1426 fprintf(stdout, "%s Sorted %d ions into separate groups of %d anions and %d cations.\n",
1427 SwS, g->nat, nAnions, nCations);
1431 /* Now we have the correct lists of anions and cations.
1432 * Copy it to the right groups.
1434 copyIndicesToGroup(indAnions, nAnions, g, cr);
1436 copyIndicesToGroup(indCations, nCations, g, cr);
1442 /*! \brief Returns TRUE if we started from an old .tpr
1444 * Then we need to re-sort anions and cations into separate groups */
1445 static gmx_bool bConvertFromOldTpr(t_swapcoords* sc)
1447 // If the last group has no atoms it means we need to convert!
1448 return (sc->ngrp >= 5) && (0 == sc->grp[4].nat);
1452 t_swap* init_swapcoords(FILE* fplog,
1453 const t_inputrec* ir,
1456 const t_state* globalState,
1457 ObservablesHistory* oh,
1459 gmx::LocalAtomSetManager* atomSets,
1460 const gmx_output_env_t* oenv,
1461 const gmx::MdrunOptions& mdrunOptions,
1462 const gmx::StartingBehavior startingBehavior)
1465 swapstateIons_t* gs;
1466 swaphistory_t* swapstate = nullptr;
1468 if ((PAR(cr)) && !DOMAINDECOMP(cr))
1470 gmx_fatal(FARGS, "Position swapping is only implemented for domain decomposition!");
1474 auto s = new t_swap();
1476 if (mdrunOptions.rerun)
1481 "%s This module does not support reruns in parallel\nPlease request a serial "
1482 "run with -nt 1 / -np 1\n",
1486 fprintf(stderr, "%s Rerun - using every available frame\n", SwS);
1488 sc->nAverage = 1; /* averaging makes no sense for reruns */
1491 if (MASTER(cr) && startingBehavior == gmx::StartingBehavior::NewSimulation)
1493 fprintf(fplog, "\nInitializing ion/water position exchanges\n");
1494 please_cite(fplog, "Kutzner2011b");
1497 switch (ir->eSwapCoords)
1499 case eswapX: s->swapdim = XX; break;
1500 case eswapY: s->swapdim = YY; break;
1501 case eswapZ: s->swapdim = ZZ; break;
1502 default: s->swapdim = -1; break;
1505 const gmx_bool bVerbose = mdrunOptions.verbose;
1507 // For compatibility with old .tpr files
1508 if (bConvertFromOldTpr(sc))
1510 convertOldToNewGroupFormat(sc, mtop, bVerbose && MASTER(cr), cr);
1513 /* Copy some data and pointers to the group structures for convenience */
1514 /* Number of atoms in the group */
1516 for (int i = 0; i < s->ngrp; i++)
1518 s->group.emplace_back(atomSets->add(
1519 gmx::ArrayRef<const int>(sc->grp[i].ind, sc->grp[i].ind + sc->grp[i].nat)));
1520 s->group[i].molname = sc->grp[i].molname;
1523 /* Check for overlapping atoms */
1524 check_swap_groups(s, mtop->natoms, bVerbose && MASTER(cr));
1526 /* Allocate space for the collective arrays for all groups */
1527 /* For the collective position array */
1528 for (int i = 0; i < s->ngrp; i++)
1531 snew(g->xc, g->atomset.numAtomsGlobal());
1533 /* For the split groups (the channels) we need some extra memory to
1534 * be able to make the molecules whole even if they span more than
1535 * half of the box size. */
1536 if ((i == eGrpSplit0) || (i == eGrpSplit1))
1538 snew(g->xc_shifts, g->atomset.numAtomsGlobal());
1539 snew(g->xc_eshifts, g->atomset.numAtomsGlobal());
1540 snew(g->xc_old, g->atomset.numAtomsGlobal());
1546 if (oh->swapHistory == nullptr)
1548 oh->swapHistory = std::make_unique<swaphistory_t>(swaphistory_t{});
1550 swapstate = oh->swapHistory.get();
1552 init_swapstate(swapstate, sc, s, mtop, globalState->x.rvec_array(), globalState->box, ir);
1555 /* After init_swapstate we have a set of (old) whole positions for our
1556 * channels. Now transfer that to all nodes */
1559 for (int ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1561 g = &(s->group[ig]);
1562 gmx_bcast((g->atomset.numAtomsGlobal()) * sizeof((g->xc_old)[0]), g->xc_old, (cr));
1566 /* Make sure that all molecules in the solvent and ion groups contain the
1567 * same number of atoms each */
1568 for (int ig = eGrpSolvent; ig < s->ngrp; ig++)
1572 g = &(s->group[ig]);
1573 g->apm = get_group_apm_check(ig, s, MASTER(cr) && bVerbose, mtop);
1575 /* Since all molecules of a group are equal, we only need enough space
1576 * to determine properties of a single molecule at at time */
1577 snew(g->m, g->apm); /* For the center of mass */
1578 charge = 0; /* To determine the charge imbalance */
1580 for (int j = 0; j < g->apm; j++)
1582 const t_atom& atom = mtopGetAtomParameters(mtop, g->atomset.globalIndex()[j], &molb);
1586 /* Total charge of one molecule of this group: */
1591 /* Need mass-weighted center of split group? */
1592 for (int j = eGrpSplit0; j <= eGrpSplit1; j++)
1595 if (sc->massw_split[j])
1597 /* Save the split group masses if mass-weighting is requested */
1598 snew(g->m, g->atomset.numAtomsGlobal());
1600 for (size_t i = 0; i < g->atomset.numAtomsGlobal(); i++)
1602 g->m[i] = mtopGetAtomMass(mtop, g->atomset.globalIndex()[i], &molb);
1607 /* Make a t_pbc struct on all nodes so that the molecules
1608 * chosen for an exchange can be made whole. */
1611 bool restartWithAppending = (startingBehavior == gmx::StartingBehavior::RestartWithAppending);
1616 fprintf(stderr, "%s Opening output file %s%s\n", SwS, fn,
1617 restartWithAppending ? " for appending" : "");
1620 s->fpout = gmx_fio_fopen(fn, restartWithAppending ? "a" : "w");
1622 if (!restartWithAppending)
1624 xvgr_header(s->fpout, "Molecule counts", "Time (ps)", "counts", exvggtXNY, oenv);
1626 for (int ig = 0; ig < s->ngrp; ig++)
1628 g = &(s->group[ig]);
1629 fprintf(s->fpout, "# %s group '%s' contains %d atom%s",
1630 ig < eSwapFixedGrpNR ? eSwapFixedGrp_names[ig] : "Ion", g->molname,
1631 static_cast<int>(g->atomset.numAtomsGlobal()),
1632 (g->atomset.numAtomsGlobal() > 1) ? "s" : "");
1633 if (!(eGrpSplit0 == ig || eGrpSplit1 == ig))
1635 fprintf(s->fpout, " with %d atom%s in each molecule of charge %g", g->apm,
1636 (g->apm > 1) ? "s" : "", g->q);
1638 fprintf(s->fpout, ".\n");
1641 fprintf(s->fpout, "#\n# Initial positions of split groups:\n");
1644 for (int j = eGrpSplit0; j <= eGrpSplit1; j++)
1647 for (size_t i = 0; i < g->atomset.numAtomsGlobal(); i++)
1649 copy_rvec(globalState->x[sc->grp[j].ind[i]], g->xc[i]);
1651 /* xc has the correct PBC representation for the two channels, so we do
1652 * not need to correct for that */
1653 get_center(g->xc, g->m, g->atomset.numAtomsGlobal(), g->center);
1654 if (!restartWithAppending)
1656 fprintf(s->fpout, "# %s group %s-center %5f nm\n", eSwapFixedGrp_names[j],
1657 DimStr[s->swapdim], g->center[s->swapdim]);
1661 if (!restartWithAppending)
1663 if ((0 != sc->bulkOffset[eCompA]) || (0 != sc->bulkOffset[eCompB]))
1665 fprintf(s->fpout, "#\n");
1667 "# You provided an offset for the position of the bulk layer(s).\n");
1669 "# That means the layers to/from which ions and water molecules are "
1672 "# are not midway (= at 0.0) between the compartment-defining layers (at "
1674 fprintf(s->fpout, "# bulk-offsetA = %g\n", sc->bulkOffset[eCompA]);
1675 fprintf(s->fpout, "# bulk-offsetB = %g\n", sc->bulkOffset[eCompB]);
1678 fprintf(s->fpout, "#\n");
1679 fprintf(s->fpout, "# Split0 cylinder radius %f nm, up %f nm, down %f nm\n", sc->cyl0r,
1680 sc->cyl0u, sc->cyl0l);
1681 fprintf(s->fpout, "# Split1 cylinder radius %f nm, up %f nm, down %f nm\n", sc->cyl1r,
1682 sc->cyl1u, sc->cyl1l);
1684 fprintf(s->fpout, "#\n");
1685 if (!mdrunOptions.rerun)
1688 "# Coupling constant (number of swap attempt steps to average over): %d "
1689 "(translates to %f ps).\n",
1690 sc->nAverage, sc->nAverage * sc->nstswap * ir->delta_t);
1691 fprintf(s->fpout, "# Threshold is %f\n", sc->threshold);
1692 fprintf(s->fpout, "#\n");
1694 "# Remarks about which atoms passed which channel use global atoms numbers "
1695 "starting at one.\n");
1704 /* Allocate memory to remember the past particle counts for time averaging */
1705 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1707 g = &(s->group[ig]);
1708 for (int ic = 0; ic < eCompNR; ic++)
1710 snew(g->comp[ic].nMolPast, sc->nAverage);
1714 /* Get the initial particle concentrations and let the other nodes know */
1717 if (startingBehavior != gmx::StartingBehavior::NewSimulation)
1719 get_initial_ioncounts_from_cpt(ir, s, swapstate, cr, bVerbose);
1723 fprintf(stderr, "%s Determining initial numbers of ions per compartment.\n", SwS);
1724 get_initial_ioncounts(ir, s, globalState->x.rvec_array(), globalState->box, cr,
1725 mdrunOptions.rerun);
1728 /* Prepare (further) checkpoint writes ... */
1729 if (startingBehavior != gmx::StartingBehavior::NewSimulation)
1731 /* Consistency check */
1732 if (swapstate->nAverage != sc->nAverage)
1734 gmx_fatal(FARGS, "%s Ion count averaging steps mismatch! checkpoint: %d, tpr: %d",
1735 SwS, swapstate->nAverage, sc->nAverage);
1740 swapstate->nAverage = sc->nAverage;
1742 fprintf(stderr, "%s Setting pointers for checkpoint writing\n", SwS);
1743 for (int ic = 0; ic < eCompNR; ic++)
1745 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1748 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1750 gs->nMolReq_p[ic] = &(g->comp[ic].nMolReq);
1751 gs->nMolPast_p[ic] = &(g->comp[ic].nMolPast[0]);
1752 gs->inflow_net_p[ic] = &(g->comp[ic].inflow_net);
1756 /* Determine the total charge imbalance */
1757 s->deltaQ = getRequestedChargeImbalance(s);
1761 fprintf(stderr, "%s Requested charge imbalance is Q(A) - Q(B) = %g e.\n", SwS, s->deltaQ);
1763 if (!restartWithAppending)
1765 fprintf(s->fpout, "# Requested charge imbalance is Q(A)-Q(B) = %g e.\n", s->deltaQ);
1771 bc_initial_concentrations(cr, ir->swap, s);
1774 /* Update the time-averaged number of molecules for all groups and compartments */
1775 for (int ig = eSwapFixedGrpNR; ig < sc->ngrp; ig++)
1778 for (int ic = 0; ic < eCompNR; ic++)
1780 update_time_window(&g->comp[ic], sc->nAverage, -1);
1784 /* Initialize arrays that keep track of through which channel the ions go */
1785 detect_flux_per_channel_init(s, swapstate, startingBehavior != gmx::StartingBehavior::NewSimulation);
1787 /* We need to print the legend if we open this file for the first time. */
1788 if (MASTER(cr) && !restartWithAppending)
1790 print_ionlist_legend(ir, s, oenv);
1796 void finish_swapcoords(t_swap* s)
1804 // Close the swap output file
1805 gmx_fio_fclose(s->fpout);
1809 /*! \brief Do we need to swap a molecule in any of the ion groups with a water molecule at this step?
1811 * From the requested and average molecule counts we determine whether a swap is needed
1812 * at this time step.
1814 static gmx_bool need_swap(t_swapcoords* sc, t_swap* s)
1819 for (ig = eSwapFixedGrpNR; ig < sc->ngrp; ig++)
1823 for (ic = 0; ic < eCompNR; ic++)
1825 if (g->comp[ic].nMolReq - g->comp[ic].nMolAv >= sc->threshold)
1835 /*! \brief Return the index of an atom or molecule suitable for swapping.
1837 * Returns the index of an atom that is far off the compartment boundaries,
1838 * that is near to the bulk layer to/from which the swaps take place.
1839 * Other atoms of the molecule (if any) will directly follow the returned index.
1841 * \param[in] comp Structure containing compartment-specific data.
1842 * \param[in] molname Name of the molecule.
1844 * \returns Index of the first atom of the molecule chosen for a position exchange.
1846 static int get_index_of_distant_atom(t_compartment* comp, const char molname[])
1849 real d = GMX_REAL_MAX;
1852 /* comp->nat contains the original number of atoms in this compartment
1853 * prior to doing any swaps. Some of these atoms may already have been
1854 * swapped out, but then they are marked with a distance of GMX_REAL_MAX
1856 for (int iMol = 0; iMol < comp->nMolBefore; iMol++)
1858 if (comp->dist[iMol] < d)
1861 d = comp->dist[ibest];
1868 "Could not get index of %s atom. Compartment contains %d %s molecules before "
1870 molname, comp->nMolBefore, molname);
1873 /* Set the distance of this index to infinity such that it won't get selected again in
1876 comp->dist[ibest] = GMX_REAL_MAX;
1878 return comp->ind[ibest];
1882 /*! \brief Swaps centers of mass and makes molecule whole if broken */
1883 static void translate_positions(rvec* x, int apm, rvec old_com, rvec new_com, t_pbc* pbc)
1886 rvec reference, dx, correctPBCimage;
1889 /* Use the first atom as the reference for PBC */
1890 copy_rvec(x[0], reference);
1892 for (i = 0; i < apm; i++)
1894 /* PBC distance between position and reference */
1895 pbc_dx(pbc, x[i], reference, dx);
1897 /* Add PBC distance to reference */
1898 rvec_add(reference, dx, correctPBCimage);
1900 /* Subtract old_com from correct image and add new_com */
1901 rvec_dec(correctPBCimage, old_com);
1902 rvec_inc(correctPBCimage, new_com);
1904 copy_rvec(correctPBCimage, x[i]);
1909 /*! \brief Write back the modified local positions from the collective array to the official positions. */
1910 static void apply_modified_positions(swap_group* g, rvec x[])
1912 auto collectiveIndex = g->atomset.collectiveIndex().begin();
1913 for (const auto localIndex : g->atomset.localIndex())
1915 /* Copy the possibly modified position */
1916 copy_rvec(g->xc[*collectiveIndex], x[localIndex]);
1922 gmx_bool do_swapcoords(t_commrec* cr,
1927 gmx_wallcycle* wcycle,
1934 int j, ic, ig, nswaps;
1935 int thisC, otherC; /* Index into this compartment and the other one */
1936 gmx_bool bSwap = FALSE;
1937 t_swapgrp * g, *gsol;
1939 rvec com_solvent, com_particle; /* solvent and swap molecule's center of mass */
1942 wallcycle_start(wcycle, ewcSWAP);
1946 set_pbc(s->pbc, ir->ePBC, box);
1948 /* Assemble the positions of the split groups, i.e. the channels.
1949 * Here we also pass a shifts array to communicate_group_positions(), so that it can make
1950 * the molecules whole even in cases where they span more than half of the box in
1952 for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1954 g = &(s->group[ig]);
1955 communicate_group_positions(cr, g->xc, g->xc_shifts, g->xc_eshifts, TRUE, x,
1956 g->atomset.numAtomsGlobal(), g->atomset.numAtomsLocal(),
1957 g->atomset.localIndex().data(),
1958 g->atomset.collectiveIndex().data(), g->xc_old, box);
1960 get_center(g->xc, g->m, g->atomset.numAtomsGlobal(), g->center); /* center of split groups == channels */
1963 /* Assemble the positions of the ions (ig = 3, 4, ...). These molecules should
1964 * be small and we can always make them whole with a simple distance check.
1965 * Therefore we pass NULL as third argument. */
1966 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1968 g = &(s->group[ig]);
1969 communicate_group_positions(cr, g->xc, nullptr, nullptr, FALSE, x, g->atomset.numAtomsGlobal(),
1970 g->atomset.numAtomsLocal(), g->atomset.localIndex().data(),
1971 g->atomset.collectiveIndex().data(), nullptr, nullptr);
1973 /* Determine how many ions of this type each compartment contains */
1974 sortMoleculesIntoCompartments(g, cr, sc, s, box, step, s->fpout, bRerun, FALSE);
1977 /* Output how many ions are in the compartments */
1980 print_ionlist(s, t, "");
1983 /* If we are doing a rerun, we are finished here, since we cannot perform
1990 /* Do we have to perform a swap? */
1991 bSwap = need_swap(sc, s);
1994 /* Since we here know that we have to perform ion/water position exchanges,
1995 * we now assemble the solvent positions */
1996 g = &(s->group[eGrpSolvent]);
1997 communicate_group_positions(cr, g->xc, nullptr, nullptr, FALSE, x, g->atomset.numAtomsGlobal(),
1998 g->atomset.numAtomsLocal(), g->atomset.localIndex().data(),
1999 g->atomset.collectiveIndex().data(), nullptr, nullptr);
2001 /* Determine how many molecules of solvent each compartment contains */
2002 sortMoleculesIntoCompartments(g, cr, sc, s, box, step, s->fpout, bRerun, TRUE);
2004 /* Save number of solvent molecules per compartment prior to any swaps */
2005 g->comp[eCompA].nMolBefore = g->comp[eCompA].nMol;
2006 g->comp[eCompB].nMolBefore = g->comp[eCompB].nMol;
2008 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
2010 g = &(s->group[ig]);
2012 for (ic = 0; ic < eCompNR; ic++)
2014 /* Determine in which compartment ions are missing and where they are too many */
2015 g->vacancy[ic] = g->comp[ic].nMolReq - g->comp[ic].nMolAv;
2017 /* Save number of ions per compartment prior to swaps */
2018 g->comp[ic].nMolBefore = g->comp[ic].nMol;
2022 /* Now actually perform the particle exchanges, one swap group after another */
2023 gsol = &s->group[eGrpSolvent];
2024 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
2028 for (thisC = 0; thisC < eCompNR; thisC++)
2030 /* Index to the other compartment */
2031 otherC = (thisC + 1) % eCompNR;
2033 while (g->vacancy[thisC] >= sc->threshold)
2035 /* Swap in an ion */
2037 /* Get the xc-index of the first atom of a solvent molecule of this compartment */
2038 isol = get_index_of_distant_atom(&gsol->comp[thisC], gsol->molname);
2040 /* Get the xc-index of a particle from the other compartment */
2041 iion = get_index_of_distant_atom(&g->comp[otherC], g->molname);
2043 get_molecule_center(&gsol->xc[isol], gsol->apm, gsol->m, com_solvent, s->pbc);
2044 get_molecule_center(&g->xc[iion], g->apm, g->m, com_particle, s->pbc);
2046 /* Subtract solvent molecule's center of mass and add swap particle's center of mass */
2047 translate_positions(&gsol->xc[isol], gsol->apm, com_solvent, com_particle, s->pbc);
2048 /* Similarly for the swap particle, subtract com_particle and add com_solvent */
2049 translate_positions(&g->xc[iion], g->apm, com_particle, com_solvent, s->pbc);
2051 /* Keep track of the changes */
2052 g->vacancy[thisC]--;
2053 g->vacancy[otherC]++;
2054 g->comp[thisC].nMol++;
2055 g->comp[otherC].nMol--;
2056 g->comp[thisC].inflow_net++;
2057 g->comp[otherC].inflow_net--;
2058 /* Correct the past time window to still get the right averages from now on */
2059 g->comp[thisC].nMolAv++;
2060 g->comp[otherC].nMolAv--;
2061 for (j = 0; j < sc->nAverage; j++)
2063 g->comp[thisC].nMolPast[j]++;
2064 g->comp[otherC].nMolPast[j]--;
2066 /* Clear ion history */
2069 int iMol = iion / g->apm;
2070 g->channel_label[iMol] = eChHistPassedNone;
2071 g->comp_from[iMol] = eDomainNotset;
2073 /* That was the swap */
2078 if (nswaps && bVerbose)
2080 fprintf(stderr, "%s Performed %d swap%s in step %" PRId64 " for iontype %s.\n", SwS,
2081 nswaps, nswaps > 1 ? "s" : "", step, g->molname);
2085 if (s->fpout != nullptr)
2087 print_ionlist(s, t, " # after swap");
2090 /* For the solvent and user-defined swap groups, each rank writes back its
2091 * (possibly modified) local positions to the official position array. */
2092 for (ig = eGrpSolvent; ig < s->ngrp; ig++)
2095 apply_modified_positions(g, x);
2098 } /* end of if(bSwap) */
2100 wallcycle_stop(wcycle, ewcSWAP);