2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017, 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"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/fileio/confio.h"
54 #include "gromacs/fileio/gmxfio.h"
55 #include "gromacs/fileio/xvgr.h"
56 #include "gromacs/gmxlib/network.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/groupcoord.h"
59 #include "gromacs/mdlib/mdrun.h"
60 #include "gromacs/mdlib/sim_util.h"
61 #include "gromacs/mdtypes/commrec.h"
62 #include "gromacs/mdtypes/inputrec.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/mdtypes/observableshistory.h"
65 #include "gromacs/mdtypes/swaphistory.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/topology/mtop_lookup.h"
69 #include "gromacs/topology/topology.h"
70 #include "gromacs/utility/cstringutil.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/pleasecite.h"
73 #include "gromacs/utility/smalloc.h"
74 #include "gromacs/utility/snprintf.h"
76 static const char *SwS = {"SWAP:"}; /**< For output that comes from the swap module */
77 static const char *SwSEmpty = {" "}; /**< Placeholder for multi-line output */
78 static const char* CompStr[eCompNR] = {"A", "B" }; /**< Compartment name */
79 static const char *SwapStr[eSwapTypesNR+1] = { "", "X-", "Y-", "Z-", nullptr}; /**< Name for the swap types. */
80 static const char *DimStr[DIM+1] = { "X", "Y", "Z", nullptr}; /**< Name for the swap dimension. */
82 /** Keep track of through which channel the ions have passed */
83 enum eChannelHistory {
84 eChHistPassedNone, eChHistPassedCh0, eChHistPassedCh1, eChHistNr
86 static const char* ChannelString[eChHistNr] = { "none", "channel0", "channel1" }; /**< Name for the channels */
88 /*! \brief Domain identifier.
90 * Keeps track of from which compartment the ions came before passing the
94 eDomainNotset, eDomainA, eDomainB, eDomainNr
96 static const char* DomainString[eDomainNr] = { "not_assigned", "Domain_A", "Domain_B" }; /**< Name for the domains */
101 * Structure containing compartment-specific data.
103 typedef struct swap_compartment
105 int nMol; /**< Number of ion or water molecules detected
106 in this compartment. */
107 int nMolBefore; /**< Number of molecules before swapping. */
108 int nMolReq; /**< Requested number of molecules in compartment. */
109 real nMolAv; /**< Time-averaged number of molecules matching
110 the compartment conditions. */
111 int *nMolPast; /**< Past molecule counts for time-averaging. */
112 int *ind; /**< Indices to collective array of atoms. */
113 real *dist; /**< Distance of atom to bulk layer, which is
114 normally the center layer of the compartment */
115 int nalloc; /**< Allocation size for ind array. */
116 int inflow_net; /**< Net inflow of ions into this compartment. */
121 * This structure contains data needed for the groups involved in swapping:
122 * split group 0, split group 1, solvent group, ion groups.
124 typedef struct swap_group
126 char *molname; /**< Name of the group or ion type */
127 int nat; /**< Number of atoms in the group */
128 int apm; /**< Number of atoms in each molecule */
129 int *ind; /**< Global atom indices of the group (size nat) */
130 int *ind_loc; /**< Local atom indices of the group */
131 int nat_loc; /**< Number of local group atoms */
132 int nalloc_loc; /**< Allocation size for ind_loc */
133 rvec *xc; /**< Collective array of group atom positions (size nat) */
134 ivec *xc_shifts; /**< Current (collective) shifts (size nat) */
135 ivec *xc_eshifts; /**< Extra shifts since last DD step (size nat) */
136 rvec *xc_old; /**< Old (collective) positions (size nat) */
137 real q; /**< Total charge of one molecule of this group */
138 int *c_ind_loc; /**< Position of local atoms in the
139 collective array, [0..nat_loc] */
140 real *m; /**< Masses (can be omitted, size apm) */
141 unsigned char *comp_from; /**< (Collective) Stores from which compartment this
142 molecule has come. This way we keep track of
143 through which channel an ion permeates
144 (size nMol = nat/apm) */
145 unsigned char *comp_now; /**< In which compartment this ion is now (size nMol) */
146 unsigned char *channel_label; /**< Which channel was passed at last by this ion?
148 rvec center; /**< Center of the group; COM if masses are used */
149 t_compartment comp[eCompNR]; /**< Distribution of particles of this group across
150 the two compartments */
151 real vacancy[eCompNR]; /**< How many molecules need to be swapped in? */
152 int fluxfromAtoB[eChanNR]; /**< Net flux of ions per channel */
153 int nCyl[eChanNR]; /**< Number of ions residing in a channel */
154 int nCylBoth; /**< Ions assigned to cyl0 and cyl1. Not good. */
159 * Main (private) data structure for the position swapping protocol.
161 typedef struct t_swap
163 int swapdim; /**< One of XX, YY, ZZ */
164 t_pbc *pbc; /**< Needed to make molecules whole. */
165 FILE *fpout; /**< Output file. */
166 int ngrp; /**< Number of t_swapgrp groups */
167 t_swapgrp *group; /**< Separate groups for channels, solvent, ions */
168 int fluxleak; /**< Flux not going through any of the channels. */
169 real deltaQ; /**< The charge imbalance between the compartments. */
174 /*! \brief Check whether point is in channel.
176 * A channel is a cylinder defined by a disc
177 * with radius r around its center c. The thickness of the cylinder is
184 * <---------c--------->
190 * \param[in] point The position (xyz) under consideration.
191 * \param[in] center The center of the cylinder.
192 * \param[in] d_up The upper extension of the cylinder.
193 * \param[in] d_down The lower extension.
194 * \param[in] r_cyl2 Cylinder radius squared.
195 * \param[in] pbc Structure with info about periodic boundary conditions.
196 * \param[in] normal The membrane normal direction is typically 3, i.e. z, but can be x or y also.
198 * \returns Whether the point is inside the defined cylindric channel.
200 static gmx_bool is_in_channel(
210 int plane1, plane2; /* Directions tangential to membrane */
213 plane1 = (normal + 1) % 3; /* typically 0, i.e. XX */
214 plane2 = (normal + 2) % 3; /* typically 1, i.e. YY */
216 /* Get the distance vector dr between the point and the center of the cylinder */
217 pbc_dx(pbc, point, center, dr); /* This puts center in the origin */
219 /* Check vertical direction */
220 if ( (dr[normal] > d_up) || (dr[normal] < -d_down) )
225 /* Check radial direction */
226 if ( (dr[plane1]*dr[plane1] + dr[plane2]*dr[plane2]) > r_cyl2)
231 /* All check passed, this point is in the cylinder */
236 /*! \brief Prints output to CompEL output file.
238 * Prints to swap output file how many ions are in each compartment,
239 * where the centers of the split groups are, and how many ions of each type
240 * passed the channels.
242 static void print_ionlist(
245 const char comment[])
248 fprintf(s->fpout, "%12.5e", time);
250 // Output number of molecules and difference to reference counts for each
251 // compartment and ion type
252 for (int iComp = 0; iComp < eCompNR; iComp++)
254 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
256 t_compartment *comp = &s->group[ig].comp[iComp];
258 fprintf(s->fpout, "%10d%10.1f%10d", comp->nMol, comp->nMolAv - comp->nMolReq, comp->inflow_net);
262 // Output center of split groups
263 fprintf(s->fpout, "%10g%10g",
264 s->group[eGrpSplit0].center[s->swapdim],
265 s->group[eGrpSplit1].center[s->swapdim]);
267 // Output ion flux for each channel and ion type
268 for (int iChan = 0; iChan < eChanNR; iChan++)
270 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
272 t_swapgrp *g = &s->group[ig];
273 fprintf(s->fpout, "%10d", g->fluxfromAtoB[iChan]);
277 /* Output the number of molecules that leaked from A to B */
278 fprintf(s->fpout, "%10d", s->fluxleak);
280 fprintf(s->fpout, "%s\n", comment);
284 /*! \brief Get the center of a group of nat atoms.
286 * Since with PBC an atom group might not be whole, use the first atom as the
287 * reference atom and determine the center with respect to this reference.
289 static void get_molecule_center(
297 rvec weightedPBCimage;
299 rvec reference, correctPBCimage, dx;
302 /* Use the first atom as the reference and put other atoms near that one */
303 /* This does not work for large molecules that span > half of the box! */
304 copy_rvec(x[0], reference);
306 /* Calculate either the weighted center or simply the center of geometry */
309 for (i = 0; i < nat; i++)
311 /* PBC distance between position and reference */
312 pbc_dx(pbc, x[i], reference, dx);
314 /* Add PBC distance to reference */
315 rvec_add(reference, dx, correctPBCimage);
317 /* Take weight into account */
318 if (nullptr == weights)
327 svmul(wi, correctPBCimage, weightedPBCimage);
330 rvec_inc(center, weightedPBCimage);
334 svmul(1.0/wsum, center, center);
339 /*! \brief Return TRUE if position x of ion (or water) is found in the compartment,
340 * i.e. between w1 and w2.
342 * One can define and additional offset "b" if one wants to exchange ions/water
343 * to or from a plane not directly in the middle of w1 and w2. The offset can be
344 * in ]-1.0, ..., +1.0 [.
345 * A bulkOffset of 0.0 means 'no offset', so the swap-layer is directly in the
346 * middle between w1 and w2. Offsets -1.0 < b < 0.0 will yield swaps nearer to w1,
347 * whereas offsets 0.0 < 0 < +1.0 will yield swaps nearer to w2.
351 * ||--------------+-------------|-------------+------------------------||
352 * w1 ? ? ? ? ? ? ? ? ? ? ? w2
353 * ||--------------+-------------|----b--------+------------------------||
358 * \param[in] w1 Position of 'wall' atom 1.
359 * \param[in] w2 Position of 'wall' atom 2.
360 * \param[in] x Position of the ion or the water molecule under consideration.
361 * \param[in] l Length of the box, from || to || in the sketch.
362 * \param[in] bulkOffset Where is the bulk layer "b" to be found between w1 and w2?
363 * \param[out] distance_from_b Distance of x to the bulk layer "b".
365 * \returns TRUE if x is between w1 and w2.
367 * Also computes the distance of x to the compartment center (the layer that is
368 * normally situated in the middle of w1 and w2 that would be considered as having
369 * the bulk concentration of ions).
371 static gmx_bool compartment_contains_atom(
377 real *distance_from_b)
383 /* First set the origin in the middle of w1 and w2 */
390 /* Now choose the PBC image of x that is closest to the origin: */
401 *distance_from_b = (real)fabs(x - bulkOffset*0.5*width);
403 /* Return TRUE if we now are in area "????" */
404 if ( (x >= w1) && (x < w2) )
415 /*! \brief Updates the time-averaged number of ions in a compartment. */
416 static void update_time_window(t_compartment *comp, int values, int replace)
422 /* Put in the new value */
425 comp->nMolPast[replace] = comp->nMol;
428 /* Compute the new time-average */
430 for (i = 0; i < values; i++)
432 average += comp->nMolPast[i];
435 comp->nMolAv = average;
439 /*! \brief Add the atom with collective index ci to the atom list in compartment 'comp'.
441 * \param[in] ci Index of this ion in the collective xc array.
442 * \param[inout] comp Compartment to add this atom to.
443 * \param[in] distance Shortest distance of this atom to the bulk layer,
444 * from which ion/water pairs are selected for swapping.
446 static void add_to_list(
453 if (nr >= comp->nalloc)
455 comp->nalloc = over_alloc_dd(nr+1);
456 srenew(comp->ind, comp->nalloc);
457 srenew(comp->dist, comp->nalloc);
460 comp->dist[nr] = distance;
465 /*! \brief Determine the compartment boundaries from the channel centers. */
466 static void get_compartment_boundaries(
470 real *left, real *right)
473 real leftpos, rightpos, leftpos_orig;
478 gmx_fatal(FARGS, "No compartment %c.", c+'A');
481 pos0 = s->group[eGrpSplit0].center[s->swapdim];
482 pos1 = s->group[eGrpSplit1].center[s->swapdim];
495 /* This gets us the other compartment: */
498 leftpos_orig = leftpos;
500 rightpos = leftpos_orig + box[s->swapdim][s->swapdim];
508 /*! \brief Determine the per-channel ion flux.
510 * To determine the flux through the individual channels, we
511 * remember the compartment and channel history of each ion. An ion can be
512 * either in channel0 or channel1, or in the remaining volume of compartment
516 * +-----------------+
519 * ||||||||||0|||||||| bilayer with channel 0
524 * |||||1||||||||||||| bilayer with channel 1
527 * +-----------------+
531 static void detect_flux_per_channel(
536 unsigned char *comp_now,
537 unsigned char *comp_from,
538 unsigned char *channel_label,
548 gmx_bool in_cyl0, in_cyl1;
555 /* Check whether ion is inside any of the channels */
556 in_cyl0 = is_in_channel(atomPosition, s->group[eGrpSplit0].center, sc->cyl0u, sc->cyl0l, cyl0_r2, s->pbc, sd);
557 in_cyl1 = is_in_channel(atomPosition, s->group[eGrpSplit1].center, sc->cyl1u, sc->cyl1l, cyl1_r2, s->pbc, sd);
559 if (in_cyl0 && in_cyl1)
561 /* Ion appears to be in both channels. Something is severely wrong! */
563 *comp_now = eDomainNotset;
564 *comp_from = eDomainNotset;
565 *channel_label = eChHistPassedNone;
569 /* Ion is in channel 0 now */
570 *channel_label = eChHistPassedCh0;
571 *comp_now = eDomainNotset;
576 /* Ion is in channel 1 now */
577 *channel_label = eChHistPassedCh1;
578 *comp_now = eDomainNotset;
583 /* Ion is not in any of the channels, so it must be in domain A or B */
586 *comp_now = eDomainA;
590 *comp_now = eDomainB;
594 /* Only take action, if ion is now in domain A or B, and was before
595 * in the other domain!
597 if (eDomainNotset == *comp_from)
599 /* Maybe we can set the domain now */
600 *comp_from = *comp_now; /* Could still be eDomainNotset, though */
602 else if ( (*comp_now != eDomainNotset ) /* if in channel */
603 && (*comp_from != *comp_now) )
605 /* Obviously the ion changed its domain.
606 * Count this for the channel through which it has passed. */
607 switch (*channel_label)
609 case eChHistPassedNone:
612 fprintf(stderr, " %s Warning! Step %s, ion %d moved from %s to %s\n",
613 SwS, gmx_step_str(step, buf), iAtom, DomainString[*comp_from], DomainString[*comp_now]);
616 fprintf(stderr, ", possibly due to a swap in the original simulation.\n");
620 fprintf(stderr, "but did not pass cyl0 or cyl1 as defined in the .mdp file.\n"
621 "Do you have an ion somewhere within the membrane?\n");
622 /* Write this info to the CompEL output file: */
623 fprintf(s->fpout, " # Warning: step %s, ion %d moved from %s to %s (probably through the membrane)\n",
624 gmx_step_str(step, buf), iAtom,
625 DomainString[*comp_from], DomainString[*comp_now]);
629 case eChHistPassedCh0:
630 case eChHistPassedCh1:
631 if (*channel_label == eChHistPassedCh0)
640 if (eDomainA == *comp_from)
642 g->fluxfromAtoB[chan_nr]++;
646 g->fluxfromAtoB[chan_nr]--;
648 fprintf(fpout, "# Atom nr. %d finished passing %s.\n", iAtom, ChannelString[*channel_label]);
651 gmx_fatal(FARGS, "%s Unknown channel history entry for ion type '%s'\n",
656 /* This ion has moved to the _other_ compartment ... */
657 *comp_from = *comp_now;
658 /* ... and it did not pass any channel yet */
659 *channel_label = eChHistPassedNone;
664 /*! \brief Determines which ions or solvent molecules are in compartment A and B */
665 static void sortMoleculesIntoCompartments(
675 gmx_swapcoords_t s = sc->si_priv;
676 int nMolNotInComp[eCompNR]; /* consistency check */
677 real cyl0_r2 = sc->cyl0r * sc->cyl0r;
678 real cyl1_r2 = sc->cyl1r * sc->cyl1r;
680 /* Get us a counter that cycles in the range of [0 ... sc->nAverage[ */
681 int replace = (step/sc->nstswap) % sc->nAverage;
683 for (int comp = eCompA; comp <= eCompB; comp++)
687 /* Get lists of atoms that match criteria for this compartment */
688 get_compartment_boundaries(comp, sc->si_priv, box, &left, &right);
690 /* First clear the ion molecule lists */
691 g->comp[comp].nMol = 0;
692 nMolNotInComp[comp] = 0; /* consistency check */
694 /* Loop over the molecules and atoms of this group */
695 for (int iMol = 0, iAtom = 0; iAtom < g->nat; iAtom += g->apm, iMol++)
700 /* Is this first atom of the molecule in the compartment that we look at? */
701 if (compartment_contains_atom(left, right, g->xc[iAtom][sd], box[sd][sd], sc->bulkOffset[comp], &dist) )
703 /* Add the first atom of this molecule to the list of molecules in this compartment */
704 add_to_list(iAtom, &g->comp[comp], dist);
706 /* Master also checks for ion groups through which channel each ion has passed */
707 if (MASTER(cr) && (g->comp_now != nullptr) && !bIsSolvent)
709 int globalAtomNr = g->ind[iAtom] + 1; /* PDB index starts at 1 ... */
710 detect_flux_per_channel(g, globalAtomNr, comp, g->xc[iAtom],
711 &g->comp_now[iMol], &g->comp_from[iMol], &g->channel_label[iMol],
712 sc, cyl0_r2, cyl1_r2, step, bRerun, fpout);
717 nMolNotInComp[comp]++;
720 /* Correct the time-averaged number of ions in the compartment */
723 update_time_window(&g->comp[comp], sc->nAverage, replace);
727 /* Flux detection warnings */
728 if (MASTER(cr) && !bIsSolvent)
733 "%s Warning: %d atoms were detected as being in both channels! Probably your split\n"
734 "%s cylinder is way too large, or one compartment has collapsed (step %" GMX_PRId64 ")\n",
735 SwS, g->nCylBoth, SwS, step);
737 fprintf(s->fpout, "Warning: %d atoms were assigned to both channels!\n", g->nCylBoth);
743 if (bIsSolvent && nullptr != fpout)
745 fprintf(fpout, "# Solv. molecules in comp.%s: %d comp.%s: %d\n",
746 CompStr[eCompA], g->comp[eCompA].nMol,
747 CompStr[eCompB], g->comp[eCompB].nMol);
750 /* Consistency checks */
751 if (nMolNotInComp[eCompA] + nMolNotInComp[eCompB] != g->nat / g->apm)
753 fprintf(stderr, "%s Warning: Inconsistency while assigning '%s' molecules to compartments. !inA: %d, !inB: %d, total molecules %d\n",
754 SwS, g->molname, nMolNotInComp[eCompA], nMolNotInComp[eCompB], g->nat/g->apm);
757 int sum = g->comp[eCompA].nMol + g->comp[eCompB].nMol;
758 if (sum != g->nat/g->apm)
760 fprintf(stderr, "%s Warning: %d molecules are in group '%s', but altogether %d have been assigned to the compartments.\n",
761 SwS, g->nat/g->apm, g->molname, sum);
766 /*! \brief Find out how many group atoms are in the compartments initially */
767 static void get_initial_ioncounts(
769 rvec x[], /* the initial positions */
785 /* Loop over the user-defined (ion) groups */
786 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
790 /* Copy the initial positions of the atoms in the group
791 * to the collective array so that we can compartmentalize */
792 for (i = 0; i < g->nat; i++)
795 copy_rvec(x[ind], g->xc[i]);
798 /* Set up the compartments and get lists of atoms in each compartment */
799 sortMoleculesIntoCompartments(g, cr, sc, box, 0, s->fpout, bRerun, FALSE);
801 /* Set initial molecule counts if requested (as signaled by "-1" value) */
802 for (ic = 0; ic < eCompNR; ic++)
804 int requested = sc->grp[ig].nmolReq[ic];
807 g->comp[ic].nMolReq = g->comp[ic].nMol;
811 g->comp[ic].nMolReq = requested;
815 /* Check whether the number of requested molecules adds up to the total number */
816 req = g->comp[eCompA].nMolReq + g->comp[eCompB].nMolReq;
817 tot = g->comp[eCompA].nMol + g->comp[eCompB].nMol;
821 gmx_fatal(FARGS, "Mismatch of the number of %s ions summed over both compartments.\n"
822 "You requested a total of %d ions (%d in A and %d in B),\n"
823 "but there are a total of %d ions of this type in the system.\n",
824 g->molname, req, g->comp[eCompA].nMolReq,
825 g->comp[eCompB].nMolReq, tot);
828 /* Initialize time-averaging:
829 * Write initial concentrations to all time bins to start with */
830 for (ic = 0; ic < eCompNR; ic++)
832 g->comp[ic].nMolAv = g->comp[ic].nMol;
833 for (i = 0; i < sc->nAverage; i++)
835 g->comp[ic].nMolPast[i] = g->comp[ic].nMol;
842 /*! \brief Copy history of ion counts from checkpoint file.
844 * When called, the checkpoint file has already been read in. Here we copy
845 * over the values from .cpt file to the swap data structure.
847 static void get_initial_ioncounts_from_cpt(
848 t_inputrec *ir, swaphistory_t *swapstate,
849 t_commrec *cr, gmx_bool bVerbose)
861 /* Copy the past values from the checkpoint values that have been read in already */
864 fprintf(stderr, "%s Copying values from checkpoint\n", SwS);
867 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
870 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
872 for (int ic = 0; ic < eCompNR; ic++)
874 g->comp[ic].nMolReq = gs->nMolReq[ic];
875 g->comp[ic].inflow_net = gs->inflow_net[ic];
879 fprintf(stderr, "%s ... Influx netto: %d Requested: %d Past values: ", SwS,
880 g->comp[ic].inflow_net, g->comp[ic].nMolReq);
883 for (int j = 0; j < sc->nAverage; j++)
885 g->comp[ic].nMolPast[j] = gs->nMolPast[ic][j];
888 fprintf(stderr, "%d ", g->comp[ic].nMolPast[j]);
893 fprintf(stderr, "\n");
901 /*! \brief The master lets all others know about the initial ion counts. */
902 static void bc_initial_concentrations(
913 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
917 for (ic = 0; ic < eCompNR; ic++)
919 gmx_bcast(sizeof(g->comp[ic].nMolReq), &(g->comp[ic].nMolReq), cr);
920 gmx_bcast(sizeof(g->comp[ic].nMol ), &(g->comp[ic].nMol ), cr);
921 gmx_bcast( swap->nAverage * sizeof(g->comp[ic].nMolPast[0]), g->comp[ic].nMolPast, cr);
927 /*! \brief Ensure that each atom belongs to at most one of the swap groups. */
928 static void check_swap_groups(t_swap *s, int nat, gmx_bool bVerbose)
930 int *nGroup = nullptr; /* This array counts for each atom in the MD system to
931 how many swap groups it belongs (should be 0 or 1!) */
933 int nMultiple = 0; /* Number of atoms belonging to multiple groups */
938 fprintf(stderr, "%s Making sure each atom belongs to at most one of the swap groups.\n", SwS);
941 /* Add one to the group count of atoms belonging to a swap group: */
943 for (int i = 0; i < s->ngrp; i++)
945 t_swapgrp *g = &s->group[i];
946 for (int j = 0; j < g->nat; j++)
948 /* Get the global index of this atom of this group: */
953 /* Make sure each atom belongs to at most one of the groups: */
954 for (int i = 0; i < nat; i++)
965 gmx_fatal(FARGS, "%s Cannot perform swapping since %d atom%s allocated to more than one swap index group.\n"
966 "%s Each atom must be allocated to at most one of the split groups, the swap groups, or the solvent.\n"
967 "%s Check the .mdp file settings regarding the swap index groups or the index groups themselves.\n",
968 SwS, nMultiple, (1 == nMultiple) ? " is" : "s are", SwSEmpty, SwSEmpty);
973 /*! \brief Get the number of atoms per molecule for this group.
975 * Also ensure that all the molecules in this group have this number of atoms.
977 static int get_group_apm_check(
983 t_swapgrp *g = &s->group[igroup];
984 int *ind = s->group[igroup].ind;
985 int nat = s->group[igroup].nat;
987 /* Determine the number of solvent atoms per solvent molecule from the
988 * first solvent atom: */
990 mtopGetMolblockIndex(mtop, ind[0], &molb, nullptr, nullptr);
991 int apm = mtop->molblock[molb].natoms_mol;
995 fprintf(stderr, "%s Checking whether all %s molecules consist of %d atom%s\n", SwS,
996 g->molname, apm, apm > 1 ? "s" : "");
999 /* Check whether this is also true for all other solvent atoms */
1000 for (int i = 1; i < nat; i++)
1002 mtopGetMolblockIndex(mtop, ind[i], &molb, nullptr, nullptr);
1003 if (apm != mtop->molblock[molb].natoms_mol)
1005 gmx_fatal(FARGS, "Not all molecules of swap group %d consist of %d atoms.",
1010 //TODO: check whether charges and masses of each molecule are identical!
1015 /*! \brief Print the legend to the swap output file.
1017 * Also print the initial values of ion counts and position of split groups.
1019 static void print_ionlist_legend(t_inputrec *ir,
1020 const gmx_output_env_t *oenv)
1022 const char **legend;
1026 t_swap *s = ir->swap->si_priv;
1027 int nIonTypes = ir->swap->ngrp - eSwapFixedGrpNR;
1028 snew(legend, eCompNR*nIonTypes*3 + 2 + eChanNR*nIonTypes + 1);
1030 // Number of molecules and difference to reference counts for each
1031 // compartment and ion type
1032 for (int ic = count = 0; ic < eCompNR; ic++)
1034 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1036 t_swapGroup *g = &ir->swap->grp[ig];
1037 real q = s->group[ig].q;
1039 snprintf(buf, STRLEN, "%s %s ions (charge %s%g)", CompStr[ic], g->molname, q > 0 ? "+" : "", q);
1040 legend[count++] = gmx_strdup(buf);
1042 snprintf(buf, STRLEN, "%s av. mismatch to %d %s ions",
1043 CompStr[ic], s->group[ig].comp[ic].nMolReq, g->molname);
1044 legend[count++] = gmx_strdup(buf);
1046 snprintf(buf, STRLEN, "%s net %s ion influx", CompStr[ic], g->molname);
1047 legend[count++] = gmx_strdup(buf);
1051 // Center of split groups
1052 snprintf(buf, STRLEN, "%scenter of %s of split group 0", SwapStr[ir->eSwapCoords], (nullptr != s->group[eGrpSplit0].m) ? "mass" : "geometry");
1053 legend[count++] = gmx_strdup(buf);
1054 snprintf(buf, STRLEN, "%scenter of %s of split group 1", SwapStr[ir->eSwapCoords], (nullptr != s->group[eGrpSplit1].m) ? "mass" : "geometry");
1055 legend[count++] = gmx_strdup(buf);
1057 // Ion flux for each channel and ion type
1058 for (int ic = 0; ic < eChanNR; ic++)
1060 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1062 t_swapGroup *g = &ir->swap->grp[ig];
1063 snprintf(buf, STRLEN, "A->ch%d->B %s permeations", ic, g->molname);
1064 legend[count++] = gmx_strdup(buf);
1068 // Number of molecules that leaked from A to B
1069 snprintf(buf, STRLEN, "leakage");
1070 legend[count++] = gmx_strdup(buf);
1072 xvgr_legend(s->fpout, count, legend, oenv);
1074 fprintf(s->fpout, "# Instantaneous ion counts and time-averaged differences to requested numbers\n");
1076 // We add a simple text legend helping to identify the columns with xvgr legend strings
1077 fprintf(s->fpout, "# time (ps)");
1078 for (int i = 0; i < count; i++)
1080 snprintf(buf, STRLEN, "s%d", i);
1081 fprintf(s->fpout, "%10s", buf);
1083 fprintf(s->fpout, "\n");
1088 /*! \brief Initialize channel ion flux detection routine.
1090 * Initialize arrays that keep track of where the ions come from and where
1093 static void detect_flux_per_channel_init(
1095 swaphistory_t *swapstate,
1096 gmx_bool bStartFromCpt)
1099 swapstateIons_t *gs;
1101 /* All these flux detection routines run on the master only */
1102 if (swapstate == nullptr)
1107 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1110 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1112 /******************************************************/
1113 /* Channel and domain history for the individual ions */
1114 /******************************************************/
1115 if (bStartFromCpt) /* set the pointers right */
1117 g->comp_from = gs->comp_from;
1118 g->channel_label = gs->channel_label;
1120 else /* allocate memory for molecule counts */
1122 snew(g->comp_from, g->nat/g->apm);
1123 gs->comp_from = g->comp_from;
1124 snew(g->channel_label, g->nat/g->apm);
1125 gs->channel_label = g->channel_label;
1127 snew(g->comp_now, g->nat/g->apm);
1129 /* Initialize the channel and domain history counters */
1130 for (int i = 0; i < g->nat/g->apm; i++)
1132 g->comp_now[i] = eDomainNotset;
1135 g->comp_from[i] = eDomainNotset;
1136 g->channel_label[i] = eChHistPassedNone;
1140 /************************************/
1141 /* Channel fluxes for both channels */
1142 /************************************/
1143 g->nCyl[eChan0] = 0;
1144 g->nCyl[eChan1] = 0;
1150 fprintf(stderr, "%s Copying channel fluxes from checkpoint file data\n", SwS);
1154 // Loop over ion types (and both channels)
1155 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1158 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1160 for (int ic = 0; ic < eChanNR; ic++)
1162 fprintf(stderr, "%s Channel %d flux history for ion type %s (charge %g): ", SwS, ic, g->molname, g->q);
1165 g->fluxfromAtoB[ic] = gs->fluxfromAtoB[ic];
1169 g->fluxfromAtoB[ic] = 0;
1172 fprintf(stderr, "%d molecule%s",
1173 g->fluxfromAtoB[ic], g->fluxfromAtoB[ic] == 1 ? "" : "s");
1174 fprintf(stderr, "\n");
1178 /* Set pointers for checkpoint writing */
1179 swapstate->fluxleak_p = &s->fluxleak;
1180 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1183 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1185 for (int ic = 0; ic < eChanNR; ic++)
1187 gs->fluxfromAtoB_p[ic] = &g->fluxfromAtoB[ic];
1193 /*! \brief Outputs the initial structure to PDB file for debugging reasons.
1195 * Output the starting structure so that in case of multimeric channels
1196 * the user can check whether we have the correct PBC image for all atoms.
1197 * If this is not correct, the ion counts per channel will be very likely
1200 static void outputStartStructureIfWanted(gmx_mtop_t *mtop, rvec *x, int ePBC, matrix box)
1202 char *env = getenv("GMX_COMPELDUMP");
1206 fprintf(stderr, "\n%s Found env.var. GMX_COMPELDUMP, will output CompEL starting structure made whole.\n"
1207 "%s In case of multimeric channels, please check whether they have the correct PBC representation.\n",
1210 write_sto_conf_mtop("CompELAssumedWholeConfiguration.pdb", *mtop->name, mtop, x, nullptr, ePBC, box);
1215 /*! \brief Initialize the swapstate structure, used for checkpoint writing.
1217 * The swapstate struct stores the information we need to make the channels
1218 * whole again after restarts from a checkpoint file. Here we do the following:
1219 * a) If we did not start from .cpt, we prepare the struct for proper .cpt writing,
1220 * b) if we did start from .cpt, we copy over the last whole structures from .cpt,
1221 * c) in any case, for subsequent checkpoint writing, we set the pointers in
1222 * swapstate to the x_old arrays, which contain the correct PBC representation of
1223 * multimeric channels at the last time step.
1225 static void init_swapstate(
1226 swaphistory_t *swapstate,
1229 rvec x[], /* the initial positions */
1233 rvec *x_pbc = nullptr; /* positions of the whole MD system with molecules made whole */
1240 /* We always need the last whole positions such that
1241 * in the next time step we can make the channels whole again in PBC */
1242 if (swapstate->bFromCpt)
1244 /* Copy the last whole positions of each channel from .cpt */
1245 g = &(s->group[eGrpSplit0]);
1246 for (int i = 0; i < g->nat; i++)
1248 copy_rvec(swapstate->xc_old_whole[eChan0][i], g->xc_old[i]);
1250 g = &(s->group[eGrpSplit1]);
1251 for (int i = 0; i < g->nat; i++)
1253 copy_rvec(swapstate->xc_old_whole[eChan1][i], g->xc_old[i]);
1258 swapstate->eSwapCoords = ir->eSwapCoords;
1260 /* Set the number of ion types and allocate memory for checkpointing */
1261 swapstate->nIonTypes = s->ngrp - eSwapFixedGrpNR;
1262 snew(swapstate->ionType, swapstate->nIonTypes);
1264 /* Store the total number of ions of each type in the swapstateIons
1265 * structure that is accessible during checkpoint writing */
1266 for (int ii = 0; ii < swapstate->nIonTypes; ii++)
1268 swapstateIons_t *gs = &swapstate->ionType[ii];
1269 gs->nMol = sc->grp[ii + eSwapFixedGrpNR].nat;
1272 /* Extract the initial split group positions. */
1274 /* Remove pbc, make molecule whole. */
1275 snew(x_pbc, mtop->natoms);
1276 copy_rvecn(x, x_pbc, 0, mtop->natoms);
1278 /* This can only make individual molecules whole, not multimers */
1279 do_pbc_mtop(nullptr, ir->ePBC, box, mtop, x_pbc);
1281 /* Output the starting structure? */
1282 outputStartStructureIfWanted(mtop, x_pbc, ir->ePBC, box);
1284 /* If this is the first run (i.e. no checkpoint present) we assume
1285 * that the starting positions give us the correct PBC representation */
1286 for (int ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1288 g = &(s->group[ig]);
1289 for (int i = 0; i < g->nat; i++)
1291 copy_rvec(x_pbc[g->ind[i]], g->xc_old[i]);
1296 /* Prepare swapstate arrays for later checkpoint writing */
1297 swapstate->nat[eChan0] = s->group[eGrpSplit0].nat;
1298 swapstate->nat[eChan1] = s->group[eGrpSplit1].nat;
1301 /* For subsequent checkpoint writing, set the swapstate pointers to the xc_old
1302 * arrays that get updated at every swapping step */
1303 swapstate->xc_old_whole_p[eChan0] = &s->group[eGrpSplit0].xc_old;
1304 swapstate->xc_old_whole_p[eChan1] = &s->group[eGrpSplit1].xc_old;
1307 /*! \brief Determine the total charge imbalance resulting from the swap groups */
1308 static real getRequestedChargeImbalance(t_swap *s)
1313 real particle_charge;
1314 real particle_number[eCompNR];
1316 // s->deltaQ = ( (-1) * s->comp[eCompA][eIonNEG].nat_req + s->comp[eCompA][eIonPOS].nat_req )
1317 // - ( (-1) * s->comp[eCompB][eIonNEG].nat_req + s->comp[eCompB][eIonPOS].nat_req );
1319 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1323 particle_charge = g->q;
1324 particle_number[eCompA] = g->comp[eCompA].nMolReq;
1325 particle_number[eCompB] = g->comp[eCompB].nMolReq;
1327 DeltaQ += particle_charge * (particle_number[eCompA] - particle_number[eCompB]);
1334 /*! \brief Sorts anions and cations into two separate groups
1336 * This routine should be called for the 'anions' and 'cations' group,
1337 * of which the indices were lumped together in the older version of the code.
1339 static void copyIndicesToGroup(
1347 /* If explicit ion counts were requested in the .mdp file
1348 * (by setting positive values for the number of ions),
1349 * we can make an additional consistency check here */
1350 if ( (g->nmolReq[eCompA] < 0) && (g->nmolReq[eCompB] < 0) )
1352 if (g->nat != (g->nmolReq[eCompA] + g->nmolReq[eCompB]) )
1354 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
1355 "%s Inconsistency while importing swap-related data from an old input file version.\n"
1356 "%s The requested ion counts in compartments A (%d) and B (%d)\n"
1357 "%s do not add up to the number of ions (%d) of this type for the group '%s'.\n",
1358 SwS, SwSEmpty, g->nmolReq[eCompA], g->nmolReq[eCompB], SwSEmpty, g->nat, g->molname);
1362 srenew(g->ind, g->nat);
1363 for (int i = 0; i < g->nat; i++)
1365 g->ind[i] = indIons[i];
1370 /*! \brief Converts old .tpr file CompEL contents to new data layout.
1372 * If we have read an old .tpr file (tpxv <= tpxv_CompElPolyatomicIonsAndMultipleIonTypes),
1373 * anions and cations are stored together in group #3. In the new
1374 * format we store each ion type in a separate group.
1375 * The 'classic' groups are:
1376 * #0 split group 0 - OK
1377 * #1 split group 1 - OK
1379 * #3 anions - contains also cations, needs to be converted
1380 * #4 cations - empty before conversion
1383 static void convertOldToNewGroupFormat(
1389 t_swapGroup *g = &sc->grp[3];
1391 /* Loop through the atom indices of group #3 (anions) and put all indices
1392 * that belong to cations into the cation group.
1396 int *indAnions = nullptr;
1397 int *indCations = nullptr;
1398 snew(indAnions, g->nat);
1399 snew(indCations, g->nat);
1402 for (int i = 0; i < g->nat; i++)
1404 const t_atom &atom = mtopGetAtomParameters(mtop, g->ind[i], &molb);
1407 // This is an anion, add it to the list of anions
1408 indAnions[nAnions++] = g->ind[i];
1412 // This is a cation, add it to the list of cations
1413 indCations[nCations++] = g->ind[i];
1419 fprintf(stdout, "%s Sorted %d ions into separate groups of %d anions and %d cations.\n",
1420 SwS, g->nat, nAnions, nCations);
1424 /* Now we have the correct lists of anions and cations.
1425 * Copy it to the right groups.
1427 copyIndicesToGroup(indAnions, nAnions, g, cr);
1429 copyIndicesToGroup(indCations, nCations, g, cr);
1437 /*! \brief Returns TRUE if we started from an old .tpr
1439 * Then we need to re-sort anions and cations into separate groups */
1440 static gmx_bool bConvertFromOldTpr(t_swapcoords *sc)
1442 // If the last group has no atoms it means we need to convert!
1443 if ( (sc->ngrp >= 5) && (0 == sc->grp[4].nat) )
1451 void init_swapcoords(
1459 ObservablesHistory *oh,
1461 const gmx_output_env_t *oenv,
1462 unsigned long Flags)
1467 swapstateIons_t *gs;
1468 gmx_bool bAppend, bStartFromCpt, bRerun;
1470 swaphistory_t *swapstate = nullptr;
1473 if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
1475 gmx_fatal(FARGS, "Position swapping is only implemented for domain decomposition!");
1478 bAppend = Flags & MD_APPENDFILES;
1479 bStartFromCpt = Flags & MD_STARTFROMCPT;
1480 bRerun = Flags & MD_RERUN;
1483 snew(sc->si_priv, 1);
1490 gmx_fatal(FARGS, "%s This module does not support reruns in parallel\nPlease request a serial run with -nt 1 / -np 1\n", SwS);
1493 fprintf(stderr, "%s Rerun - using every available frame\n", SwS);
1495 sc->nAverage = 1; /* averaging makes no sense for reruns */
1498 if (MASTER(cr) && !bAppend)
1500 fprintf(fplog, "\nInitializing ion/water position exchanges\n");
1501 please_cite(fplog, "Kutzner2011b");
1504 switch (ir->eSwapCoords)
1520 // For compatibility with old .tpr files
1521 if (bConvertFromOldTpr(sc) )
1523 convertOldToNewGroupFormat(sc, mtop, bVerbose && MASTER(cr), cr);
1526 /* Copy some data and pointers to the group structures for convenience */
1527 /* Number of atoms in the group */
1529 snew(s->group, s->ngrp);
1530 for (int i = 0; i < s->ngrp; i++)
1532 s->group[i].nat = sc->grp[i].nat;
1533 s->group[i].ind = sc->grp[i].ind;
1534 s->group[i].molname = sc->grp[i].molname;
1537 /* Check for overlapping atoms */
1538 check_swap_groups(s, mtop->natoms, bVerbose && MASTER(cr));
1540 /* Allocate space for the collective arrays for all groups */
1541 /* For the collective position array */
1542 for (int i = 0; i < s->ngrp; i++)
1545 snew(g->xc, g->nat);
1546 snew(g->c_ind_loc, g->nat);
1548 /* For the split groups (the channels) we need some extra memory to
1549 * be able to make the molecules whole even if they span more than
1550 * half of the box size. */
1551 if ( (i == eGrpSplit0) || (i == eGrpSplit1) )
1553 snew(g->xc_shifts, g->nat);
1554 snew(g->xc_eshifts, g->nat);
1555 snew(g->xc_old, g->nat);
1561 if (oh->swapHistory == nullptr)
1563 oh->swapHistory = std::unique_ptr<swaphistory_t>(new swaphistory_t {});
1565 swapstate = oh->swapHistory.get();
1567 init_swapstate(swapstate, sc, mtop, x, box, ir);
1570 /* After init_swapstate we have a set of (old) whole positions for our
1571 * channels. Now transfer that to all nodes */
1574 for (int ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1576 g = &(s->group[ig]);
1577 gmx_bcast((g->nat)*sizeof((g->xc_old)[0]), g->xc_old, (cr));
1581 /* Make sure that all molecules in the solvent and ion groups contain the
1582 * same number of atoms each */
1583 for (int ig = eGrpSolvent; ig < s->ngrp; ig++)
1587 g = &(s->group[ig]);
1588 g->apm = get_group_apm_check(ig, s, MASTER(cr) && bVerbose, mtop);
1590 /* Since all molecules of a group are equal, we only need enough space
1591 * to determine properties of a single molecule at at time */
1592 snew(g->m, g->apm); /* For the center of mass */
1593 charge = 0; /* To determine the charge imbalance */
1595 for (int j = 0; j < g->apm; j++)
1597 const t_atom &atom = mtopGetAtomParameters(mtop, g->ind[j], &molb);
1601 /* Total charge of one molecule of this group: */
1606 /* Need mass-weighted center of split group? */
1607 for (int j = eGrpSplit0; j <= eGrpSplit1; j++)
1610 if (TRUE == sc->massw_split[j])
1612 /* Save the split group masses if mass-weighting is requested */
1615 for (int i = 0; i < g->nat; i++)
1617 g->m[i] = mtopGetAtomMass(mtop, g->ind[i], &molb);
1622 /* Make a t_pbc struct on all nodes so that the molecules
1623 * chosen for an exchange can be made whole. */
1625 /* Every node needs to call set_pbc() and therefore every node needs
1626 * to know the box dimensions */
1627 copy_mat(box, boxCopy);
1630 gmx_bcast(sizeof(boxCopy), boxCopy, cr);
1632 set_pbc(s->pbc, -1, boxCopy);
1638 fprintf(stderr, "%s Opening output file %s%s\n", SwS, fn, bAppend ? " for appending" : "");
1641 s->fpout = gmx_fio_fopen(fn, bAppend ? "a" : "w" );
1645 xvgr_header(s->fpout, "Molecule counts", "Time (ps)", "counts", exvggtXNY, oenv);
1647 for (int ig = 0; ig < s->ngrp; ig++)
1649 g = &(s->group[ig]);
1650 fprintf(s->fpout, "# %s group '%s' contains %d atom%s",
1651 ig < eSwapFixedGrpNR ? eSwapFixedGrp_names[ig] : "Ion",
1652 g->molname, g->nat, (g->nat > 1) ? "s" : "");
1653 if (!(eGrpSplit0 == ig || eGrpSplit1 == ig) )
1655 fprintf(s->fpout, " with %d atom%s in each molecule of charge %g",
1656 g->apm, (g->apm > 1) ? "s" : "", g->q);
1658 fprintf(s->fpout, ".\n");
1661 fprintf(s->fpout, "#\n# Initial positions of split groups:\n");
1664 for (int j = eGrpSplit0; j <= eGrpSplit1; j++)
1667 for (int i = 0; i < g->nat; i++)
1669 copy_rvec(x[sc->grp[j].ind[i]], g->xc[i]);
1671 /* xc has the correct PBC representation for the two channels, so we do
1672 * not need to correct for that */
1673 get_center(g->xc, g->m, g->nat, g->center);
1676 fprintf(s->fpout, "# %s group %s-center %5f nm\n", eSwapFixedGrp_names[j],
1677 DimStr[s->swapdim], g->center[s->swapdim]);
1683 if ( (0 != sc->bulkOffset[eCompA]) || (0 != sc->bulkOffset[eCompB]) )
1685 fprintf(s->fpout, "#\n");
1686 fprintf(s->fpout, "# You provided an offset for the position of the bulk layer(s).\n");
1687 fprintf(s->fpout, "# That means the layers to/from which ions and water molecules are swapped\n");
1688 fprintf(s->fpout, "# are not midway (= at 0.0) between the compartment-defining layers (at +/- 1.0).\n");
1689 fprintf(s->fpout, "# bulk-offsetA = %g\n", sc->bulkOffset[eCompA]);
1690 fprintf(s->fpout, "# bulk-offsetB = %g\n", sc->bulkOffset[eCompB]);
1693 fprintf(s->fpout, "#\n");
1694 fprintf(s->fpout, "# Split0 cylinder radius %f nm, up %f nm, down %f nm\n",
1695 sc->cyl0r, sc->cyl0u, sc->cyl0l);
1696 fprintf(s->fpout, "# Split1 cylinder radius %f nm, up %f nm, down %f nm\n",
1697 sc->cyl1r, sc->cyl1u, sc->cyl1l);
1699 fprintf(s->fpout, "#\n");
1702 fprintf(s->fpout, "# Coupling constant (number of swap attempt steps to average over): %d (translates to %f ps).\n",
1703 sc->nAverage, sc->nAverage*sc->nstswap*ir->delta_t);
1704 fprintf(s->fpout, "# Threshold is %f\n", sc->threshold);
1705 fprintf(s->fpout, "#\n");
1706 fprintf(s->fpout, "# Remarks about which atoms passed which channel use global atoms numbers starting at one.\n");
1715 /* Prepare for parallel or serial run */
1718 for (int ig = 0; ig < s->ngrp; ig++)
1720 g = &(s->group[ig]);
1723 g->ind_loc = nullptr;
1728 for (int ig = 0; ig < s->ngrp; ig++)
1730 g = &(s->group[ig]);
1731 g->nat_loc = g->nat;
1732 g->ind_loc = g->ind;
1733 /* c_ind_loc needs to be set to identity in the serial case */
1734 for (int i = 0; i < g->nat; i++)
1736 g->c_ind_loc[i] = i;
1741 /* Allocate memory to remember the past particle counts for time averaging */
1742 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1744 g = &(s->group[ig]);
1745 for (int ic = 0; ic < eCompNR; ic++)
1747 snew(g->comp[ic].nMolPast, sc->nAverage);
1751 /* Get the initial particle concentrations and let the other nodes know */
1756 get_initial_ioncounts_from_cpt(ir, swapstate, cr, bVerbose);
1760 fprintf(stderr, "%s Determining initial numbers of ions per compartment.\n", SwS);
1761 get_initial_ioncounts(ir, x, box, cr, bRerun);
1764 /* Prepare (further) checkpoint writes ... */
1767 /* Consistency check */
1768 if (swapstate->nAverage != sc->nAverage)
1770 gmx_fatal(FARGS, "%s Ion count averaging steps mismatch! checkpoint: %d, tpr: %d",
1771 SwS, swapstate->nAverage, sc->nAverage);
1776 swapstate->nAverage = sc->nAverage;
1778 fprintf(stderr, "%s Setting pointers for checkpoint writing\n", SwS);
1779 for (int ic = 0; ic < eCompNR; ic++)
1781 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1784 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1786 gs->nMolReq_p[ic] = &(g->comp[ic].nMolReq);
1787 gs->nMolPast_p[ic] = &(g->comp[ic].nMolPast[0]);
1788 gs->inflow_net_p[ic] = &(g->comp[ic].inflow_net);
1792 /* Determine the total charge imbalance */
1793 s->deltaQ = getRequestedChargeImbalance(s);
1797 fprintf(stderr, "%s Requested charge imbalance is Q(A) - Q(B) = %g e.\n", SwS, s->deltaQ);
1801 fprintf(s->fpout, "# Requested charge imbalance is Q(A)-Q(B) = %g e.\n", s->deltaQ);
1807 bc_initial_concentrations(cr, ir->swap);
1810 /* Update the time-averaged number of molecules for all groups and compartments */
1811 for (int ig = eSwapFixedGrpNR; ig < sc->ngrp; ig++)
1814 for (int ic = 0; ic < eCompNR; ic++)
1816 update_time_window(&g->comp[ic], sc->nAverage, -1);
1820 /* Initialize arrays that keep track of through which channel the ions go */
1821 detect_flux_per_channel_init(s, swapstate, bStartFromCpt);
1823 /* We need to print the legend if we open this file for the first time. */
1824 if (MASTER(cr) && !bAppend)
1826 print_ionlist_legend(ir, oenv);
1831 void finish_swapcoords(t_swapcoords *sc)
1833 if (sc->si_priv->fpout)
1835 // Close the swap output file
1836 gmx_fio_fclose(sc->si_priv->fpout);
1841 void dd_make_local_swap_groups(gmx_domdec_t *dd, t_swapcoords *sc)
1847 /* Make split groups, solvent group, and user-defined groups of particles
1848 * under control of the swap protocol */
1849 for (ig = 0; ig < sc->ngrp; ig++)
1851 g = &(sc->si_priv->group[ig]);
1852 dd_make_local_group_indices(dd->ga2la, g->nat, g->ind,
1853 &(g->nat_loc), &(g->ind_loc), &(g->nalloc_loc), g->c_ind_loc);
1858 /*! \brief Do we need to swap a molecule in any of the ion groups with a water molecule at this step?
1860 * From the requested and average molecule counts we determine whether a swap is needed
1861 * at this time step.
1863 static gmx_bool need_swap(t_swapcoords *sc)
1871 for (ig = eSwapFixedGrpNR; ig < sc->ngrp; ig++)
1875 for (ic = 0; ic < eCompNR; ic++)
1877 if (g->comp[ic].nMolReq - g->comp[ic].nMolAv >= sc->threshold)
1887 /*! \brief Return the index of an atom or molecule suitable for swapping.
1889 * Returns the index of an atom that is far off the compartment boundaries,
1890 * that is near to the bulk layer to/from which the swaps take place.
1891 * Other atoms of the molecule (if any) will directly follow the returned index.
1893 * \param[in] comp Structure containing compartment-specific data.
1894 * \param[in] molname Name of the molecule.
1896 * \returns Index of the first atom of the molecule chosen for a position exchange.
1898 static int get_index_of_distant_atom(
1899 t_compartment *comp,
1900 const char molname[])
1903 real d = GMX_REAL_MAX;
1906 /* comp->nat contains the original number of atoms in this compartment
1907 * prior to doing any swaps. Some of these atoms may already have been
1908 * swapped out, but then they are marked with a distance of GMX_REAL_MAX
1910 for (int iMol = 0; iMol < comp->nMolBefore; iMol++)
1912 if (comp->dist[iMol] < d)
1915 d = comp->dist[ibest];
1921 gmx_fatal(FARGS, "Could not get index of %s atom. Compartment contains %d %s molecules before swaps.",
1922 molname, comp->nMolBefore, molname);
1925 /* Set the distance of this index to infinity such that it won't get selected again in
1928 comp->dist[ibest] = GMX_REAL_MAX;
1930 return comp->ind[ibest];
1934 /*! \brief Swaps centers of mass and makes molecule whole if broken */
1935 static void translate_positions(
1943 rvec reference, dx, correctPBCimage;
1946 /* Use the first atom as the reference for PBC */
1947 copy_rvec(x[0], reference);
1949 for (i = 0; i < apm; i++)
1951 /* PBC distance between position and reference */
1952 pbc_dx(pbc, x[i], reference, dx);
1954 /* Add PBC distance to reference */
1955 rvec_add(reference, dx, correctPBCimage);
1957 /* Subtract old_com from correct image and add new_com */
1958 rvec_dec(correctPBCimage, old_com);
1959 rvec_inc(correctPBCimage, new_com);
1961 copy_rvec(correctPBCimage, x[i]);
1966 /*! \brief Write back the the modified local positions from the collective array to the official positions. */
1967 static void apply_modified_positions(
1974 for (l = 0; l < g->nat_loc; l++)
1976 /* Get the right local index to write to */
1978 /* Where is the local atom in the collective array? */
1979 cind = g->c_ind_loc[l];
1981 /* Copy the possibly modified position */
1982 copy_rvec(g->xc[cind], x[ii]);
1987 gmx_bool do_swapcoords(
1992 gmx_wallcycle *wcycle,
2000 int j, ic, ig, nswaps;
2001 int thisC, otherC; /* Index into this compartment and the other one */
2002 gmx_bool bSwap = FALSE;
2003 t_swapgrp *g, *gsol;
2005 rvec com_solvent, com_particle; /* solvent and swap molecule's center of mass */
2008 wallcycle_start(wcycle, ewcSWAP);
2014 /* Assemble the positions of the split groups, i.e. the channels.
2015 * Here we also pass a shifts array to communicate_group_positions(), so that it can make
2016 * the molecules whole even in cases where they span more than half of the box in
2018 for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
2020 g = &(s->group[ig]);
2021 communicate_group_positions(cr, g->xc, g->xc_shifts, g->xc_eshifts, TRUE,
2022 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, g->xc_old, box);
2024 get_center(g->xc, g->m, g->nat, g->center); /* center of split groups == channels */
2027 /* Assemble the positions of the ions (ig = 3, 4, ...). These molecules should
2028 * be small and we can always make them whole with a simple distance check.
2029 * Therefore we pass NULL as third argument. */
2030 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
2032 g = &(s->group[ig]);
2033 communicate_group_positions(cr, g->xc, nullptr, nullptr, FALSE,
2034 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, nullptr, nullptr);
2036 /* Determine how many ions of this type each compartment contains */
2037 sortMoleculesIntoCompartments(g, cr, sc, box, step, s->fpout, bRerun, FALSE);
2040 /* Output how many ions are in the compartments */
2043 print_ionlist(s, t, "");
2046 /* If we are doing a rerun, we are finished here, since we cannot perform
2053 /* Do we have to perform a swap? */
2054 bSwap = need_swap(sc);
2057 /* Since we here know that we have to perform ion/water position exchanges,
2058 * we now assemble the solvent positions */
2059 g = &(s->group[eGrpSolvent]);
2060 communicate_group_positions(cr, g->xc, nullptr, nullptr, FALSE,
2061 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, nullptr, nullptr);
2063 /* Determine how many molecules of solvent each compartment contains */
2064 sortMoleculesIntoCompartments(g, cr, sc, box, step, s->fpout, bRerun, TRUE);
2066 /* Save number of solvent molecules per compartment prior to any swaps */
2067 g->comp[eCompA].nMolBefore = g->comp[eCompA].nMol;
2068 g->comp[eCompB].nMolBefore = g->comp[eCompB].nMol;
2070 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
2072 g = &(s->group[ig]);
2074 for (ic = 0; ic < eCompNR; ic++)
2076 /* Determine in which compartment ions are missing and where they are too many */
2077 g->vacancy[ic] = g->comp[ic].nMolReq - g->comp[ic].nMolAv;
2079 /* Save number of ions per compartment prior to swaps */
2080 g->comp[ic].nMolBefore = g->comp[ic].nMol;
2084 /* Now actually perform the particle exchanges, one swap group after another */
2085 gsol = &s->group[eGrpSolvent];
2086 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
2090 for (thisC = 0; thisC < eCompNR; thisC++)
2092 /* Index to the other compartment */
2093 otherC = (thisC+1) % eCompNR;
2095 while (g->vacancy[thisC] >= sc->threshold)
2097 /* Swap in an ion */
2099 /* Get the xc-index of the first atom of a solvent molecule of this compartment */
2100 isol = get_index_of_distant_atom(&gsol->comp[thisC], gsol->molname);
2102 /* Get the xc-index of a particle from the other compartment */
2103 iion = get_index_of_distant_atom(&g->comp[otherC], g->molname);
2105 get_molecule_center(&gsol->xc[isol], gsol->apm, gsol->m, com_solvent, s->pbc);
2106 get_molecule_center(&g->xc[iion], g->apm, g->m, com_particle, s->pbc);
2108 /* Subtract solvent molecule's center of mass and add swap particle's center of mass */
2109 translate_positions(&gsol->xc[isol], gsol->apm, com_solvent, com_particle, s->pbc);
2110 /* Similarly for the swap particle, subtract com_particle and add com_solvent */
2111 translate_positions(&g->xc[iion], g->apm, com_particle, com_solvent, s->pbc);
2113 /* Keep track of the changes */
2114 g->vacancy[thisC ]--;
2115 g->vacancy[otherC]++;
2116 g->comp [thisC ].nMol++;
2117 g->comp [otherC].nMol--;
2118 g->comp [thisC ].inflow_net++;
2119 g->comp [otherC].inflow_net--;
2120 /* Correct the past time window to still get the right averages from now on */
2121 g->comp [thisC ].nMolAv++;
2122 g->comp [otherC].nMolAv--;
2123 for (j = 0; j < sc->nAverage; j++)
2125 g->comp[thisC ].nMolPast[j]++;
2126 g->comp[otherC].nMolPast[j]--;
2128 /* Clear ion history */
2131 int iMol = iion / g->apm;
2132 g->channel_label[iMol] = eChHistPassedNone;
2133 g->comp_from[iMol] = eDomainNotset;
2135 /* That was the swap */
2140 if (nswaps && bVerbose)
2142 fprintf(stderr, "%s Performed %d swap%s in step %" GMX_PRId64 " for iontype %s.\n",
2143 SwS, nswaps, nswaps > 1 ? "s" : "", step, g->molname);
2147 if (s->fpout != nullptr)
2149 print_ionlist(s, t, " # after swap");
2152 /* For the solvent and user-defined swap groups, each rank writes back its
2153 * (possibly modified) local positions to the official position array. */
2154 for (ig = eGrpSolvent; ig < s->ngrp; ig++)
2157 apply_modified_positions(g, x);
2160 } /* end of if(bSwap) */
2162 wallcycle_stop(wcycle, ewcSWAP);